1. Functions overview

1.1. Tracking of a particule and operations allong a flow line

class icetrackpy.icetracker.Mesh(mesh_path=None, compiled=True, parallel=True, nb_processors=8)

Bases: object

This object contains the mesh on wich we want to work and all of the values associated with it, it can also add values to the mesh by for example computing the strain field from the velocity field

This methode creats an object of Mesh, is is used at the begining to initialise the mesh object

Parameters:
  • mesh_path (str or None) – This is the complete path of the mesh on which you want to work on, it should contain the name of the file and the extension (idealy a .vtu), for example it can be something like : ‘/home/my_mesh.vtu’, if no values or given the function “from_save” can be used to restor a saved mesh

  • compiled (bool) – Defines if later the object can use compiled versions of the code or if the python versions of the algorihtms should be used

  • parallel (bool) – Defines if multiprocessing is allowed or not later when dooing operations on the mesh

  • nb_processors (int) – If parallel==True, it defines in how many parts the task is splited, this is idealy equal to the number of cores of the processor

DN(a1: float, a2: float, a3: float)

This function gives out the derivative of the shape functions, this enables the uder to then calculate the jacobian for example Theses shape functions are for wedges elements only The numeroation is made acording to this page https://web.mit.edu/calculix_v2.7/CalculiX/ccx_2.7/doc/ccx/node34.html

Parameters:
  • a1 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 1

  • a2 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 2

  • a3 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 3

Returns:

returns a matrix of size 6 by 3 of the derivative of the shape functions

Return type:

np.array([6,3])

J(a1: float, a2: float, a3: float, wedge_id: int)

This function calculates the jacobian of the wedge element at the cordinate given and for the given wedge

Parameters:
  • a1 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 1

  • a2 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 2

  • a3 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 3

  • wedge_id (int) – the id of the wedge in wich the jacobian will be calculated

Returns:

returns a matrix of size 3 by 3 which is the jacobian

Return type:

np.array([3,3])

N(a1, a2, a3)

This methode computes the shape functions of a wedge element base on the numeroation on this web site ‘https://web.mit.edu/calculix_v2.7/CalculiX/ccx_2.7/doc/ccx/node34.html

Parameters:
  • a1 (float) – coordiate in the 1 direction in the local bass

  • a2 (float) – coordiate in the 2 direction in the local bass

  • a3 (float) – coordiate in the 3 direction in the local bass

Returns:

returns an array of size 6 of the shape functions

Return type:

np.array

_compute_deviatoric_stress_cuffey(strain_rate, A1, A2, Q1, Q2, T, glen_exponent)

Compute the deviatoric stress tensor directly from the strain rate tensor using Glen’s flow law with cuffey prefactor

Parameters:
  • strain_rate – Strain rate tensor (3x3 ndarray)

  • A1 – flow rate prefactor dependant on the temeprature for low temp <-10

  • A2 – flow rate prefactor dependant on the temeprature for high temp >-10

  • Q1 – activation energie for the first A1

  • Q2 – activation energie for the first A2

  • glen_exponent – glen exponent of the flow law

Returns:

Deviatoric stress tensor (3x3 ndarray).

_compute_deviatoric_stress_glen(strain_rate, A, n)

Compute the deviatoric stress tensor directly from the strain rate tensor using Glen’s flow law.

Parameters:
  • strain_rate – Strain rate tensor (3x3 ndarray)

  • A – flow rate viscosity

  • n – glen exponent

Returns:

Deviatoric stress tensor (3x3 ndarray).

_invert_mesh_parallel(queue, task)

This methode allows us to find the wedges connected to a given lis of points This methode is internal and allows for parallel computing internally, it should not be directly used

Parameters:
  • queue (mp.queue object) – queue for multiprocessing

  • task (list) – list of the start node and end node for the mesh insertion task

_invert_mesh_parallel_compiled(queue, task, mesh_wedges_ctype)

This methode allows us to find the wedges connected to a given lis of points This methode is internal and allows for parallel computing internally with the compiled function, it should not be directly used

Parameters:
  • queue (mp.queue object) – queue for multiprocessing

  • task (list) – list of the start node and end node for the mesh insertion task

  • mesh_wedges_ctype – memory adress where the list of wedges is storded

Type:

ctype object

_strain_rate_parallel(queue, tasks)

This function in not meant to used directly by the end user but it allows for strain calculation in parallel

Parameters:
  • queue (mp.queue object) – queue for multiprocessing

  • task (list) – list of the start node and end node for the mesh insertion task

compute_strain_rate(save_to_mesh=True)

This functions is meant to compute the strain rate if only the velocity field is availabe is can save or not the new filed to the vtu file

Parameters:

save_to_mesh (bool) – tells if we whant to save the calculated value to the vtu file for future use or visualisation

compute_stress_cuffey(A1=9125353404000.0, A2=7.660165593599999e+23, Q1=60000.0, Q2=115000.0, T=273, glen_exponent=3)

This function computes an approxiamtion of the stress using the glen law it enable the user to compute stress after the simulation has run, this doent igve the same results has elmer a closer look should be taken to better understand why

Parameters:
  • A1 (float) – flow rate prefactor dependant on the temeprature for low temp <-10

  • A2 (float) – flow rate prefactor dependant on the temeprature for high temp >-10

  • Q1 (float) – activation energie for the first A1

  • Q2 (float) – activation energie for the first A2

  • T (float) – temperature of the ice in kelins

  • glen_exponent (float) – exponent of the glen law, usually 3 for ice

Returns:

adds the fileds “sxx” – “syz” to the mesh nodes fields so that they can be advected with the particule

compute_stress_glen(glen_exponent=3, viscosity=3.17e-06)

This function computes an approxiamtion of the stress using the glen law it enable the user to compute stress after the simulation has run, this doent igve the same results has elmer a closer look should be taken to better understand why

Parameters:
  • glen_exponent (float) – exponent of the glen law, usually 3 for ice

  • viscosity (float) – the viscosity that should be used for the glen law

Returns:

adds the fileds “sxx” – “syz” to the mesh nodes fields so that they can be advected with the particule

extract_data()

This function extracts all of the available data at each nodes of the mesh and adds it to the local mesh object

from_save(folder_path)

Allows the user to creat a mesh from a previous save made before, it will recreat all of the elements needed for the mesh creation

Parameters:

folder_path (str) – path to the global folder of the save

grad_v(a1: float, a2: float, a3: float, wedge_id: int)

This function calculates the gradient of the velocity field in the given element at the given coordinates

Parameters:
  • a1 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 1

  • a2 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 2

  • a3 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 3

  • wedge_id (int) – the id of the wedge in wich the jacobian will be calculated

Returns:

returns a matrix of size 3 by 3 which is the gradient of the velocity

Return type:

np.array([3,3])

invert_mesh()

This methode allows the user to find the relation between the nodes and the wedges, this is not stored in the .vtu file. It is mendatory to call it after mesh import for the tracker to work properlly This operation is quite time consuming (O(n^2) with n the number of elements) thus it is parallelized and compiled

is_top_surface(point)

This function tells the user if a given point is on the top surface of the mesh or not

Parameters:

point (int) – this should be the id of the point where we whant to know if we are on top or not

Returns:

returns true if the point given is on the top surface

Return type:

bool

select_points(methode='z_min_max', ice_mask=False, ice_mask_key='icymask', ice_mask_inside_val=1, tol=0.1, **kwargs)

This methode enable the user to to select a big number of points at the same time, the selected particules can only be selected on the top and bottom surface not on the edeges or in the volume, each time they or selected at the node

Parameters:
  • methode (str) – this parameters helps select the methode of point selection it can take 2 values : “z_min_max” or “z_min_max_y_min_max_x_min_max”

  • ice_mask (bool) – says if there is an ice mask to use for point selection

  • ice_mask_key (str) – key word in the mesh dictionnary for the icemask layer

  • ice_mask_inside_val (float) – value where ice is supposed to exist

  • tol (float) – this is by how much the ice_mask_value is able to change for it to be inside the glacier or not

  • kwargs (kwargs) – this is to define the boundaries of the selection depending on the methode, if the methode “z_min_max” is selected, the arguments are “z_min”:int, “z_max”:int and “surf_type”:str in [“top”,”bottom”,”top_and_bottom”], if the methode “z_min_max_y_min_max_x_min_max” is selected the parameters “x_min”:int, “x_max”:int, “y_min”:int, “y_max”:int are added

Returns:

it returns the list of the coodrinates of each points and the list of the trianlges of the mesh composed by those points

Return type:

tuple(list,list)

strain_rate(a1: float, a2: float, a3: float, wedge_id)

This function calculates the strain rate a given coordinate in a given wedge

Parameters:
  • a1 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 1

  • a2 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 2

  • a3 (float) – value of the coordinate of the point in the local bases of the wedge for the direction 3

  • wedge_id (int) – the id of the wedge in wich the jacobian will be calculated

Returns:

returns a matrix of size 3 by 3 which is the tensor of strain rate

Return type:

np.array([3,3])

class icetrackpy.icetracker.ParticuleTracker(mesh, eps_top=0.0001, u_min=0.01, start_acu_level=0, nb_friends_max=3, ellipse_factor=10, outside_iter=1, ini_iter=100, dt=-0.01, n=1000, compiled=True, params_to_track=[], parallel=True, verbose=True)

Bases: object

This class is the instance that enables the user to do some particule tracking and more opperations on each particules

This creats an instance of particule tracker, many parameters are set at the beegining that will be used mater for the particule tracking instance

Parameters:
  • mesh (Mesh.object) – this should be an object mesh created by icetrackpy, it can be just a mesh openned by meshio for example

  • eps_top (float) – this defines the level at which a particule died, if the projection of the nomalized vector motion on the normal to the surface is higher then eps_top then the particule died, if set high it allows the particule to do some surfing, should be between 0 an 1

  • u_min (float) – This defines the minimum speed at which a particule can go, if the particule is slower then said speed it will die, the speed is both calculated by taking the norm of the velocity field and also by calculating the real speed between 2 cordiantes

  • start_acu_level – this is the level at which the particules are allowed to melt under if time goes in the normal direction and fall if time is reversed

  • nb_friends_max (int) – this defines the number of elements the code will look thru when looking for an wedge arround a point, this size is the small axis of an ellips good values between 2 and 4, smaller is faster but less accurate

  • ellipse_factor (int) – this defines how our research radius looks like an ellipse, for meshesh where elements can have a really bad shapes it is a good idea to set the value higher, value between 1 and 10 is good higher is slower

  • outside_iter (int) – this is the number of random iterations added to the normal iteration to find the good value of wedge

  • ini_iter (int) – number of initial tryins to find the node on the first step time, if the mesh is really non concave this value need to be high arround 100 is a good start

  • dt (float) – this is the step time of the particule tracker higher is faster but less precise

  • n (int) – number of maximum iterations of the particule tracker

  • compiled (bool) – tels the programm if he is able to use compiled code recources

  • params_to_track (list[str]) – list of strings of the parameters that should be added to the dictionnary during the tracking

  • parallel (bool) – tels the code if he is allowded to multi process the task

  • verbose (bool) – tells if time steps are printed or not

_flow_path_parallel(queue, point: array, particule_id)

Called function for parallelisation of the tracking shloud no be used directly by the end user

Parameters:
  • queue (queue.object) – queue of the storted results:

  • point (np.array) – point to be computed

_int_strain_particule_parallel(queue, particule: dict, param_to_int: str)

This function calls int_strain_particule, it is used for parallization it should not be used by the end user

Parameters:
  • particule (dict) – this is the dictionnary of the particule element

  • param_to_int (str) – parameter wich we want to integrate like “e”

compact_tensor(field)

This function helps the user compact a tensor that is in the form of many components in the form of xx or _11

Parameters:

filed – filed name of the tensor that should be compacted for example “e”

Returns:

adds to the said filed the values of the said tensor

compact_vector(field)

This function helps the user compact a vector that is in the form of many components in the form of x or _1

Parameters:

filed – filed name of the vector that should be compacted for example “velocity”

Returns:

adds to the said filed the values of the said vector

compute_age()

adds the filed I1 that represent the integration of the time step over the time, this represents the age of the particule

Returns:

adds the age of the particule in the field “I1”

compute_deformation_eigval(strain_eigval: str)

This function is used when the strain has been calculated and that we know the eigen value of the strain, in this situation we can compute the déformation tensor of the particule, this is usfull to use the anisotropy factor

Parameters:

strain_eigval (str) – this is the name of the field of the eigen values of the strain, this field should be a 1D array

Returns:

adds the field deformation_eigval to the particule trunk where the eigen values of the deformation tensor are storded a 1D np array

compute_dev(field: str)

This function computes the deviatoric version of the tensor field given

Parameters:

field – type of tensor field for which we want to know the deviatoric version, for example “e”

Type:

“str”

Returns:

adds to the particule dictionnary each filed of the tensor with the addition of the letter “D” for deviatoric, for example “Dexx”

compute_flatness_anisotropy(tensor_eigval)

This function computes the flatness anisotropy of a tensor field

Parameters:

tensor_eigval (str) – eigen values field of the tensor we want to calculate

Returns:

retruns the flatness anisotropy at each step time in the field flatness_anisotropy

compute_flinn(field, top_cutoff=1)

Computes the vector of the field in the flinn plane for easy visualisation of the anisotropy

Parameters:
  • field (str) – field name taht will go through the flinn process calculation, idealy it is the field of the aigen values of the transfomation vector

  • top_cutoff (int) – to help visualisation, this parametters compresses the to value, if set to one 1 every thing has the proper value, if set to 2 the maximum bound is set to half of the highest value, value shoud be an int higher then 1

compute_fractional_anisotropy(tensor_eigval)

This function computes the fractional anisotropy of a tensor field

Parameters:

tensor_eigval (str) – eigen values field of the tensor we want to calculate

Returns:

retruns the fractional anisotropy at each step time in the field relative_anisotropy

compute_in_particule_coord(tensor_type: str)

This function function helps us convert ensy tensor field that we have on the particule flowline to the local bases

Parameters:

tensor_type (str) – give what field we should convert to the local base for strain rate it will be “e” for example

Returns:

adds the fields of each component of the tensor with a “L” before to indicate the local carrater of the field

compute_path(points: array)

This functions calculates the tracking of a list of points, it is parallelised and in part compiled

Parameters:

points (list(np.array([3]))) – list of the starting points for the tracking

compute_relative_anisotropy(tensor_eigval)

This function computes the relative anisotropy of a tensor field

Parameters:

tensor_eigval (str) – eigen values field of the tensor we want to calculate

Returns:

retruns the relative anisotropy at each step time in the field relative_anisotropy

compute_rotation_matrix()

This functions computes the roation matrix from the field “path” it takes the vector at each time step and calculate the angle relative to base coordinate system. The time step vectors are filtered to keep the algorythm stable even with very noisy data. It can be used with a reversed time step or a foward time step

compute_strain(param_to_int)

This function calculates the strain rate in the high deformation context, it is the parallelised or not version depending on the object config

Parameters:

param_to_int (str) – parameter wich we want to integrate like “e”

Returns:

adds to the diconnary 2 fileds that are “eigval_I”+param to int and “eigvect_I”+param wich are the rigen values and aigen vectors of the strain for each time step

compute_volume_ratio_anisotropy(tensor_eigval)

This function computes the volume_ratio anisotropy of a tensor field

Parameters:

tensor_eigval (str) – eigen values field of the tensor we want to calculate

Returns:

retruns the volume_ratio anisotropy at each step time in the field relative_anisotropy

diag_tensor_field(field)

This function takes field of symmetric tensor and creats 2 new keys in the dictionnary with the eigen values and the eigen vectors of the the field tensor

Parameters:

field (str) – field to diagonalize

Retrun:

add to the dictionnary of the particule 2 fields that are “eigval_”+field and “eigvect_”+field

extract_scalar(field)

This function extracts the data that are formed in tensors or vectors into all of the bases components and adds it to the field of the particule

Parameters:

field (str) – field of tensor or vector that should be extracted

Retrun:

addts to the particule dictionnary each element with a “_1” or “_11” at the end to telle the user where the value comes from

extract_vector(field, direction=1)

This function extracts the data that are formed in tensors into all of the bases vectors and adds it to the field of the particule

Parameters:
  • field (str) – field of tensor or vector that should be extracted

  • direction (int) – defines wether we extract line : 0 or if we extract line : 1

Retrun:

adds to the particule dictionnary each element with a “_1” at the end to telle the user where the value comes from

find_best_sommet(coord_objectif, start_wedge)

Given the coordinate of a point and the coordinates of a wedge this function finds the nearest node of the mesh to the coordinate given This function is subject to compilation

Parameters:
  • coord_objectif (np.array([3])) – coordinate of the traget point

  • start_wedge (int) – id of the id of the start of the search

Returns:

returns the id of the best node of the mesh

Return type:

int

find_wedge_near_point(coord_objectif: array, start_wedge=None)

This function finds the best wedge for a given point, if the point is inside the wedge what will be returned is the wedge over the point, eles it will return the nerest wedge This funtion may call compiled functions

Parameters:
  • coord_objectif (np.array([3])) – cooridnates of the point we want to know the nearest wedge

  • start_wedge (int of None) – gives a initial idea of where we are in the mesh

Returns:

the wedge id in which we are the closes are inside, if we are inside or non, the local coodrinates of the given wedge

Return type:

int, bool, np.array([3])

flow_path(particule_coord, particule_id=None)

This function cumputs the flow path of an initial particule

Parameters:

particule_coord (np.arrray([3])) – inital coordinates of the particule

Returns:

returns a dictionnary with the following fields, “status”: tels what is the status of the particule, “is_inside”: bool if the particule is inside the glacier are touching a wall, “path” : coordiantes of the particule at each step, “velocity” velocity of the particule at each time step, “time_step” : time step of the simulaiton at each step, add all of the additional fields asked when the object was created

Return type:

dic

from_save(folder_path)

This methode allows the user to restor previous particule tracking instances that have run previouslly

Parameters:

path (str) – path to the folder where every_thing is saved, something like /home/savetracking_folder

int_strain_particule(particule: dict, param_to_int: str)

This function computes the strain tensor in the high deformation context

Parameters:
  • particule (dict) – particule dictionnary of which we want to compute the strain

  • param_to_int (str) – parameter wich we want to integrate like “e”

Returns:

returns the eigen values of the strain tensor, the eigen vectors of the strain tensor at each time step and the initial start coordinate

Return type:

np.array, np.array, np.array

int_strain_rate_over_path_small_displacement(param_to_int: str)

This funciton computes the strain rate with the small displacement approximation which is not correct in many cases

Parameters:

param_to_int (str) – parameter we want to icompute like “exx”

Returns:

adds to the particule dic the parameter with a “I” before for integrated like “Iexx”

is_inside_wedge(point_of_interest, id_wedge)

This methode helps to know weather a point is inside a given wedge, it is an essantial funtion in the particule tracking code

Parameters:
  • point_of_interest (int) – point id for which we desiere to know wether it is inside the wedge

  • id_wedge (int) – id of the wedge for which we whant to know if we are inside or not

Returns:

return where it is inside or not, a list of the coordinates in the local bases and a score to know how far ethe point is fro mthe given wedge

Return type:

bool, list(floats), float

normalise_vector(field)

This function is helpfull to nomalize a filed when we want to plot it where it needs to be nomed to 1

Parameters:

field (str) – field to normalise

Returns:

adds to the particule dic the vector field with a “N” before to indicate the normalisation

roation_matrix(xyz_component: array)

This function computes the rotation matrix that is needed to go from the x axis to the normal of the flow line

Parameters:

xyz_component (np.array([3,n]) or np.array([3])) – components of the vector we want to to know the orientation

Retrun:

returns an arrray of roation matricies or simply the rotation matrix depending on the input

Return type:

np.array([n,3,3]) or np.array([3,3])

save_tracking(path, name, things_to_save='all', header=None, particule_to_save='all', save_mesh=True, auto_delete=True)

Allows the user to save the data of a particule are all of them and save the usefull parameters of said particule in one convinient folder, make shure that all of the elements of your particules are floats, np.arrays (up to 2 dimentions) or lists

Parameters:
  • path (str) – path to the folder where the file should be saved

  • name (str) – name of the folder where everything will be saved

  • things_to_save (str ("all") or list of str or None) – what field should be saved for each partticule

  • header (None or list of str) – additionnal text to put allong the parameter name for more explicit description

  • particule_to_save (int or str ("all") of None or list) – wich particule is saved

  • save_mesh (bool) – saves the mesh for easy import

  • auto_delete (bool) – delets the file if it already exists, if set to False and the file is there it will print an error

time_int_over_path(param_to_int: str)

This function does time integration over the path for a parameter that is defined allong the flow line

Parameters:

param_to_int (str) – type of parameter we whant to integrate for example “1” if we have addted a list of 1ones in the particule dic and want to kwon the time

Returns:

adds to the particule dic the parameter with a “I” before for integrated

time_inversion()

This function must be called when negative time is used to flip elements in the dictionnary to make shure that for post processing it looks like time was going forward

class icetrackpy.icetracker.Plot(cmap='hsv', bgcolor='white', min_val=None, max_val=None, projection_2D=False, canvas_size=(1200, 800), projection_type='steriographic', cwheel_shadding='classic', cwheel_symetry=True, projection_base=array([[1, 0, 0], [0, 1, 0], [0, 0, 1]]), cbar_symmetric=False)

Bases: object

This methode helps the user to plot the different meshes and particules allong there flow line with the color representing vectors or saclars

This creats the object that will be plotted it, creats the window of the plot, thus a display option needs to be present on the computure for this to work

Parameters:
  • cmap (str) – this is the color map wanted for the plotting of the colors =, only available for the color bar not for the vector plot, by default set to “hsv”

  • bgcolor (str) – color of the backgroud of the plot by default “white”

  • canvas_size (tuple) – size of the canvas, higer will give better results when saving the figure

  • min_val (float) – minimum value for the color bar by default set to none and will be calculated during plotting

  • max_val (float) – maximum value for the color bar by default set to none and will be calculated during plotting

  • projection_2D (bool) – this tells if we want to project the plot in the xy plane for 2D mapping

  • projection_type ("steriographic" or "equal_angle" or "equal_area") – allows for all of the different types of projections

  • cwheel_shadding ("classic" or "none" or "full" or "full_sqrt" or "full_sqrt") – allows for different visual of the color ball and wheel, the sqrt version is to make the dark area easyer to read

  • cwheel_symetry (bool) – allows to chose if we want a symmetric or not colorwheel

  • projection_base (np.array) – bases from which the angles are calculated for vector plot, usfull if we want to change the orientation of the cwheel

  • cbar_symmetric (bool) – tells if we want the color bar to have symmetric values arround 0, usfful when using a diverging colormap

_download_image(url, output_file)

Download an image from the given URL and save it to a file.

Parameters:
  • url – URL of the image to download.

  • output_file – Path where the image will be saved.

_on_key_press(event)

This function is internal and is called when a key is pressed, it saves the image when “ctrl+s” is pressed

Parameters:

event (all) – is usless can be anything

_on_mouse_click(event)

This function is internal is regulates when a mousse click is pressed, this is for the selection of the line

Parameters:

event (all) – useless

_photo_and_close(event)

This function is internal and closes the window

Parameters:

event (all) – is usless can be anything

angle_from_vect(vect)

This function calculates the rotations that are needed to to convert a vector to the bases represented by the value Mesh.projection_base

Parameters:

vect (np.array or lsit) – vector for which we want to know the angles

Retrun:

returns the angles of the roation in radiants

Return type:

float,float

angles_to_color(angles)

This function converts 2 angles and gives the color corresponding on the color wheel chosen

Parameters:

angles (list or array) – the 2 angles for which we want to know the color

Returns:

returs an array of the color in the rbg space

Return type:

np.array

axis(size=150, border_color=None)

This function plots the axis xyz of the main plot on the bottom let side

Parameters:
  • size (int) – size of the view box of the axis

  • border_color (str) – color of the border of the axis by default set to None

background_geotiff(dem_crs, image_path)

This function helps the user to plot a geotiff image in the background of the 2 plot of the glacier, the geotiff can be from any sources and in theory in any coordinate system, the image should allways be put in the right place

Parameters:
  • dem_crs (str) – this is the coordinate system that is used for the glacier dem (the one used in elmer) for example for argentirere example it is ‘EPSG:27592’

  • dem_crs – this is the complete path to the geotiff file, the file should be a .tif only

cball(size=200, border_color=None)

This function plots a color ball at the bottom of the frame that rotates with the motion of the main figure

Parameters:
  • size (int) – size of the view box of the colorball

  • border_color (str) – color of the border of the colorball by default set to None

cbar(title='', number_of_digits=4)

This function prints the color bar on the right side of the plot

Parameters:
  • title (str) – title name for the plot, if not defined it will be the data that is plotted

  • number_of_digits (int) – number of decimals that we keep on the displayed values

cquarter(size=150, border_color=None)

This function plots the quarter of a colorwheel at the bottom right of the plot

Parameters:
  • size (int) – size of the view box of the colorwheel

  • border_color (str) – color of the border of the colorwheel by default set to None

cwheel(size=75, border_color=None)

This function plots the colorwheel at the bottom right of the plot

Parameters:
  • size (int) – size of the view box of the colorwheel

  • border_color (str) – color of the border of the colorwheel by default set to None

find_min_max(liste_particule_data, plot_field)

This function sets the min and max value of the color bar if the value has not been defined by the user to th value minimum and maximum value of the given field

Parameters:
  • liste_particule_data (list(dict)) – list of all the particule dictionnaries

  • plot_field (str) – name of the filed we want to set the maximum and minimum value to

line_selection()

This function activates the mousse selection of a flow line, when the scroll wheel button is pressed on a flow line it will change its color and print id of the particule

mesh(mesh, mesh_color='black', mesh_width=0.1)

This function plots the mesh frame, it only plots the top and the bottom surface of the mesh for better visibility

Parameters:
  • mesh (Mesh.object) – Mesh object from icetrackpy

  • mesh_color (str) – color of the mesh line, by default it is “black”

  • mesh_width (float) – size of the mesh line by default set to 0.1

mesh_full(mesh, mesh_color=(1, 1, 1, 0.2), mesh_width=0.1)

This function plots the mesh frame, it only plots the top and the bottom surface of the mesh for better visibility it plot a background of the mesh to mask goegraphic data

Parameters:
  • mesh (Mesh.object) – Mesh object from icetrackpy

  • mesh_color (str) – color of the mesh line, by default it is “black”

  • mesh_width (float) – size of the mesh line by default set to 0.1

particule_scalar(liste_particule_data, plot_field, linewidth=10)

This function plot the particule path with the data plotted as the color of the path

Parameters:
  • liste_particule_data (list(dict)) – list of all the particule dictionnaries

  • plot_field (str) – name of the filed we want to plot allong the flow line

  • linewidth (float) – size of the lien of the particule, set to default to 10

particule_vector(liste_particule_data, plot_field, linewidth=10)

This function plots the flow line with the data of a vector

Parameters:
  • liste_particule_data (list) – list of dictionnary of the particule data to plot

  • plot_field – field that we want to plot, this should be a vector field

  • linewidth – size of the line

Type:

str

set_camera_range(mesh)

This function sets the range of the camera in function of the mesh so that every point is in frame

Parameters:

mesh (Mesh.object) – Mesh object from icetrackpy

show(path='', auto_close=False, save_fig=False, auto_delete=True, key_save=True)

This function finelizes the plot and makes the figure appear for the user

Parameters:
  • path (str) – full path to the image and where it needs to be saved and the extension like : “home/img.png”

  • show (bool) – tells if we want the view to automaticly close once the rendering is done

  • save_fig (bool) – tells if we want to save the figure or not automaticaly

  • auto_delete (bool) – this parameter tells if we want to errase the picture that has the same name when we save a new one

Param:

surface_scalar(start_points_list: list, triangle_list: list, liste_particule_data: list, plot_field='is_inside', time_step=-1)

This function plots a surface of scalars values on the surface of the mesh and interpolates between each nodes

Parameters:
  • start_points_list (list) – list of points that represent a node of the mesh on which we kown the value

  • triangle_list (list) – list of trizngles that represent the mesh on which we plot everything

  • liste_particule_data (list) – list of all the particule dictionnaries

  • plot_field (str) – field that we went to represent

  • time_step (int) – time moment at wich we want to plot the value on the mesh

surface_vector(start_points_list, triangle_list, liste_particule_data, plot_field, time_step=-1)

This function plots a surface of vectors values on the surface of the mesh and interpolates between each nodes

Parameters:
  • start_points_list (list) – list of points that represent a node of the mesh on which we kown the value

  • triangle_list (list) – list of trizngles that represent the mesh on which we plot everything

  • liste_particule_data (list) – list of all the particule dictionnaries

  • plot_field (str) – field that we went to represent this should be a vector field it should be a normalised vector of dimention 3

  • time_step (int) – time moment at wich we want to plot the value on the mesh

icetrackpy.icetracker.pvtu_to_vtu(file_path)

This function converts a pvtu file into a vtu file that can later be read by the library

Parameters:

file_path (str or list) – it can be a string of the complet path of a single file or a list of files names

Returns:

save the files with the same names as before but in the vtu format, and it retruns a list of the new files names

Return type:

list

1.2. Creating a simple example of a glacier

class icetrackpy.theoretical_glacier.Glacier(vx=1.0, vz=-1.0, mesh_factor=[1, 1, 1], square_size=[50, 10, 10], type='3d')

Bases: object

This class is able to creat a theoretical glacier where all of the values of tensors fields are mathematically known.

The functions creats an object of Glacier.

Parameters:
  • vx (float) – velocity of the particule in the x direction, (m/year)

  • vz (float) – velocity of the particule in the z direction, (m/year)

  • mesh_factor (list) – compaction of the mesh in the 3 directions in numbers of nodes for 1 meter, should be a list of size 3 of int types elements bigger or equal to 1. [Default] :code:’[1,1,1]’

  • square_size (list) – size of the 3D volume where the particules will evolve, should be a list of size 3 of floats, the volume cant be negative or null. [Default] :code:’[50,10,10]’

  • type (str) – type of glacier that will be used for creating the fields on the mesh, 2 options are available, “compression” or “3d”, the first is juste a single bar in compression and the second one is a more realistic approximation of a glacier

__dict__ = mappingproxy({'__module__': 'icetrackpy.theoretical_glacier', '__doc__': '\n    This class is able to creat a theoretical glacier where all of the values of tensors fields are mathematically known.\n\n    ', '__init__': <function Glacier.__init__>, 'v': <function Glacier.v>, 'strain_rate': <function Glacier.strain_rate>, 'phi': <function Glacier.phi>, 'strain': <function Glacier.strain>, 'vtu': <function Glacier.vtu>, 'msh': <function Glacier.msh>, '__dict__': <attribute '__dict__' of 'Glacier' objects>, '__weakref__': <attribute '__weakref__' of 'Glacier' objects>, '__annotations__': {}})
__init__(vx=1.0, vz=-1.0, mesh_factor=[1, 1, 1], square_size=[50, 10, 10], type='3d')

The functions creats an object of Glacier.

Parameters:
  • vx (float) – velocity of the particule in the x direction, (m/year)

  • vz (float) – velocity of the particule in the z direction, (m/year)

  • mesh_factor (list) – compaction of the mesh in the 3 directions in numbers of nodes for 1 meter, should be a list of size 3 of int types elements bigger or equal to 1. [Default] :code:’[1,1,1]’

  • square_size (list) – size of the 3D volume where the particules will evolve, should be a list of size 3 of floats, the volume cant be negative or null. [Default] :code:’[50,10,10]’

  • type (str) – type of glacier that will be used for creating the fields on the mesh, 2 options are available, “compression” or “3d”, the first is juste a single bar in compression and the second one is a more realistic approximation of a glacier

__module__ = 'icetrackpy.theoretical_glacier'
__weakref__

list of weak references to the object (if defined)

msh(path: str, file_name='mesh')

Creats a .msh file containing the mesh that can be converted as elmer-ice input

Parameters:
  • path (str) – path where the file needs to be saved

  • file_name (str) – name of the vtu file

Note

saves a msh file at the path given and with the given file name

phi(x0: float, y0: float, z0: float, t: float)

Gives out the position of a particule with an initial position at a given time, this function is defined everywhere even out of the domain

Parameters:
  • x0 (float) – initial position of the particule in the x direction

  • y0 (float) – initial position of the particule in the y direction

  • z0 (float) – initial position of the particule in the z direction

  • t (float) – time at which we want to know the particule position can be positive of negative

Returns:

gives out a tuple of floats of size 3 giving out the 3 coordinates of the position of the particule

Return type:

tuple

strain(x0: float, y0: float, z0: float, t: float)

Gives out the strain of a particule with an initial position at a given time, this function is defined everywhere even out of the domain

Parameters:
  • x0 (float) – initial position of the particule in the x direction

  • y0 (float) – initial position of the particule in the y direction

  • z0 (float) – initial position of the particule in the z direction

  • t (float) – time at which we want to know the particule strain can be positive of negative

Returns:

gives out a matrix of size 3 representing the symmetric strain tensor

Return type:

np.array

strain_rate(x: float, y: float, z: float)

Gives out the strain rate field at any place in the 3D space, this function is defined everywhere even out of the domain

Parameters:
  • x (float) – cordinate in the x axis to evaluate the strain rate field

  • y (float) – cordinate in the y axis to evaluate the strain rate field

  • z (float) – cordinate in the z axis to evaluate the strain rate field

Returns:

gives out a square matrix of size 3 which is the symmetric strain tensor

Return type:

np.array

v(x: float, y: float, z: float)

Gives out the velocity field at any place in the 3D space,this function is defined everywhere even out of the domain

Parameters:
  • x (float) – cordinate in the x axis to evaluate the velocity field

  • y (float) – cordinate in the y axis to evaluate the velocity field

  • z (float) – cordinate in the z axis to evaluate the velocity field

Returns:

gives out a tuple of floats of length 3 contaning the velocity filed in the x,y,z directions

Return type:

tuple

vtu(path: str, file_name='mesh')

Creats a .vtu file containing the mesh and the vlaues of strain rate and velocity at each nodes

Parameters:
  • path (str) – path where the file needs to be saved

  • file_name (str) – name of the vtu file

Note

saves a vtu file at the path given and with the given file name

1.3. Strain calculation from elementray strains

icetrackpy.diag_matrix_iterative.eig_DM(D: array, M: array)

This functions allows to find a good approximation of the eigne values and eigen vectors of the product D by M where the determinant of D and M or 1 and D is diagonal

Parameters:
  • D (np.array) – Diagonal matrix of size 3 and of determiant equal to 1

  • M (np.array) – Matrix of size 3 by 3 that should be symmetric and of determinant 1

Returns:

gives out the eigen values and eigen vectors of the product DM

Return type:

np.array,np.array

icetrackpy.diag_matrix_iterative.int_strain(strain_list: list)

Calculates the cumulated strain from a list of elementry small strains that should be summed in the incomrpessible case

Parameters:

strain_list (list) – list of strain tensor (represented by 3 by 3 arrays) idealy they are incompressible strain rates tensor multiplied by a small step time

Returns:

returns a list of eigen values and an orthogonal matrix, when taking the product PDP^T we get the cumulated strain

Return type:

np.array,np.array