vortex_detection module¶
-
class
IMTreatment.vortex_detection.vortex_detection.
CritPoints
(unit_time='s')[source]¶ Bases:
object
Class representing a set of critical point associated to a VectorField.
Parameters: unit_time (string or Unit object) – Unity for the time. -
add_point
(foc=None, foc_c=None, node_i=None, node_o=None, sadd=None, time=None)[source]¶ Add a new point to the CritPoints object.
Parameters: - foc_c, node_i, node_o, sadd (foc,) – Representing the critical points at this time.
- time (number) – Time.
-
break_trajectories
(x=None, y=None, smooth_size=5, inplace=False)[source]¶ Break the trajectories in pieces. Usefull to do before using ‘get_mean_trajectory’
Parameters: - y (x,) – Position of breaking If not specified, break the trajectories where dx/dv=0
- smooth_size (number) – Smoothing used for the dx/dv=0 detection
-
change_unit
(axe, new_unit)[source]¶ Change the unit of an axe.
Parameters: - axe (string) – ‘y’ for changing the y axis unit ‘x’ for changing the x axis unit ‘time’ for changing the time unit
- new_unit (Unum.unit object or string) – The new unit.
-
clean_traj
(min_nmb_in_traj)[source]¶ Remove some isolated points.
Parameters: min_nmb_in_traj (integer) – Trajectories that have less than this amount of points are deleted. Associated points are also deleted.
-
compute_traj
(epsilon=None, close_traj=False)[source]¶ Compute cp trajectory from cp positions.
Parameters: - epsilon (number, optional) – Maximal distance between two successive points. default value is Inf.
- close_traj (bool) – If ‘True’, try to close the trajectories (better to get the cp fusion position)
-
crop
(intervx=None, intervy=None, intervt=None, ind=False, inplace=False)[source]¶ crop the point field.
-
display
(fields=None, cpkw={}, lnkw={}, display_traj=False, buffer_size=100, **kwargs)[source]¶ Display some critical points.
Parameters: fields (TemporalVectorFields object, optional) – If specified, critical points are displayed on the given field. critical lines are also computed and displayed
-
display_animate
(TF, **kw)[source]¶ Display an interactive windows with the velocity fields from TF and the critical points.
- TF : TemporalFields
- .
- kw : dict, optional
- Additional arguments for ‘TF.display()’
-
display_arch
(indice=None, time=None, field=None, cpkw={}, lnkw={})[source]¶ Display some critical points.
Parameters: - time (number, optional) – If specified, critical points associated to this time are displayed.
- indice (integer, optional) – If specified, critical points associated to this indice are displayed.
- field (VectorField object, optional) – If specified, critical points are displayed on the given field. critical lines are also computed and displayed
-
display_traj
(data='default', filt=None, **kw)[source]¶ Display the stored trajectories.
Parameters: - data (string) – If ‘default’, trajectories are plotted in a 2-dimensional plane. If ‘x’, x position of cp are plotted against time. If ‘y’, y position of cp are plotted against time.
- filt (array of boolean) – Filter on CP types.
- kw (dict, optional) – Arguments passed to plot.
-
display_traj_len_repartition
()[source]¶ display profiles with trajectories length repartition (usefull to choose the right epsilon)
-
get_mean_trajectory
(cp_type, min_len=20, min_nmb_to_avg=10, rel_diff_epsilon=0.03, rel_len_epsilon=0.25, verbose=False)[source]¶ Return a mean trajectory (based on a set of trajectories)
Parameters: - cp_type (string in ['foc', 'foc_c', 'node_i', 'node_o', 'sadd']) – or a set of trajectories Trajectory type to average.
- min_len (number) – Ignore trajectories with length smaller.
- min_nmb_to_avg (number) – Minimum number of values necessary to make an average
- rel_diff_epsilon (number) – Maximum relative difference used to determine the different mean trajectories.
- rel_len_epsilon (number) – Maximum relative length difference used to determine the different mean trajectories.
Returns: - mean_traj (MeanTrajectory object) – .
- skipped_traj (tuple of Points objects) – Skipped trajectories
Notes
-You may want to separate your lonf trajectories with ‘break_trajectories’ before using this function -Checking trajectory resemblance ake into account length of the
trajectories and integrale of their differences.
-
get_points_density
(kind, bw_method=None, resolution=100, output_format='normalized')[source]¶ Return the presence density map for the given point type.
Parameters: - kind (string) – Type of critical point for the density map (can be ‘foc’, ‘foc_c’, ‘sadd’, ‘node_i’, ‘node_o’)
- bw_method (str, scalar or callable, optional) – The method used to calculate the estimator bandwidth. This can be ‘scott’, ‘silverman’, a scalar constant or a callable. If a scalar, this will be used as std (it should aprocimately be the size of the density node you want to see). If a callable, it should take a gaussian_kde instance as only parameter and return a scalar. If None (default), ‘scott’ is used. See Notes for more details.
- resolution (integer or 2x1 tuple of integers, optional) – Resolution for the resulting field. Can be a tuple in order to specify resolution along x and y.
- output_format (string, optional) –
‘normalized’ (default) : give position probability (integral egal 1). ‘absolute’ : sum of integral over all points density egal 1. ‘ponderated’ : give position probability ponderated by the number
or points (integral egal number of points).’concentration’ : give local concentration (in point per surface).
-
get_traj_direction_changement
(cp_type, direction, smoothing=None)[source]¶ Return the position of trajectories direction changement.
Parameters: - cp_type (string in ['foc', 'foc_c', 'node_i', 'node_o', 'sadd']) – CP type to use
- direction (string in ['x', 'y']) – Direction along which look
- smoothing (number, optional) – Smoothing size (performed on value before gradient search).
Returns: chg_pts_1, chg_pts_2 – .
Return type: Points objects
-
iter
¶
-
iter_traj
¶
-
refine_cp_position
(cp_type, fields, inplace=True, verbose=True, extrema='max')[source]¶ Refine the position of the critical points by putting them on the given scalar field extrema.
Parameters: - cp_type (string in ['foc', 'foc_c', 'sadd', 'node_i', 'node_o']) – Critical points type to refine.
- fields (TemporalScalarFields object) – fields where to search for extrema.
- extrema (string in ['min', 'max']) – If ‘max’, cp are displaced on field maxima, if ‘min’, cp are displaced on field minima.
- inplace (boolean) –
.
- verbose (boolean) –
.
-
remove_point
(time=None, indice=None)[source]¶ Remove some critical points.
Parameters: - time (number, optional) – If specified, critical points associated to this time are removed.
- indice (integer, optional) – If specified, critical points associated to this indice are removed.
-
scale
(scalex=1.0, scaley=1.0, scalev=1.0, inplace=False)[source]¶ Change the scale of the axis.
Parameters: - scaley, scalev (scalex,) – scales along x, y and v
- inplace (boolean, optional) – If ‘True’, scaling is done in place, else, a new instance is returned.
-
smooth_traj
(tos='uniform', size=None)[source]¶ Smooth the CP trajectories.
Parameters: - tos (string, optional) – Type of smoothing, can be ‘uniform’ (default) or ‘gaussian’ (See ndimage module documentation for more details)
- size (number, optional) – Size of the smoothing (is radius for ‘uniform’ and sigma for ‘gaussian’). Default is 3 for ‘uniform’ and 1 for ‘gaussian’.
-
topo_simplify
(dist_min, kind='replacement')[source]¶ Simplify the topological points field.
Parameters: - dist_min (number) – Minimal distance between two points in the simplified field.
- kind (string, optional) – Algorithm used to simplify the field, can be ‘replacement’(default) for iterative replacement with center of mass or ‘only_delete’ to eliminate only the first order associates.
Returns: simpl_CP – Simplified topological points field.
Return type: CritPoints object
-
unit_time
¶
-
unit_x
¶
-
unit_y
¶
-
-
class
IMTreatment.vortex_detection.vortex_detection.
MeanTrajectory
(xy=array([], shape=(0, 2), dtype=float64), time=[], assoc_real_times=[], assoc_std_x=[], assoc_std_y=[], nmb_traj_used=0, base_trajs=[], unit_x='', unit_y='', unit_times='', name='')[source]¶ Bases:
IMTreatment.core.points.Points
-
crop
(intervx=None, intervy=None, intervt=None, inplace=True, ind=False)[source]¶ Crop the points cloud.
Parameters: - intervx (2x1 tuple) – Interval on x axis
- intervy (2x1 tuple) – Interval on y axis
- intervt (2x1 tuple) – Interval on time
Returns: tmp_pts – croped version of the point cloud.
Return type: Points object
-
display
(*args, **kwargs)[source]¶ Display the set of points.
Parameters: - kind (string, optional) – Can be ‘plot’ (default if points have not values). or ‘scatter’ (default if points have values). or ‘colored_plot’.
- axe_y, axe_color (axe_x,) – To determine wich value has to be plotted along which axis, and whith value is used to color the scattered points. Default plot ‘y’ to ‘x’ with colors from ‘v’.
- **plotargs (dic) – Additionnal arguments sent to ‘plot’ or ‘scatter’
-
display_error_bars
(kind='bar', **kwargs)[source]¶ Display the error bar according to the trajectories used to compute the mean trajectory.
Parameters: kind (string in ['bar', 'envelope']) – If ‘bar’, display error bar, if ‘envelope’, display envelope around trajectory
-
reconstruct_fields
(TF)[source]¶ Do a conditionnal averaging based on the CP positions.
Parameters: TF (TemporalFields) – .
-
scale
(scalex=1.0, scaley=1.0, scalet=1.0, inplace=False)[source]¶ Change the scale of the axis.
Parameters: - scaley, scalev (scalex,) – scales along x, y and v
- inplace (boolean, optional) – If ‘True’, scaling is done in place, else, a new instance is returned.
-
time
¶
-
unit_times
¶
-
used_pts_number
¶
-
used_traj_number
¶
-
-
class
IMTreatment.vortex_detection.vortex_detection.
TopoPoints
[source]¶ Bases:
object
Represent a topological points field. (only for topo simplification, in fact, should be merged with CritPoints)
Parameters: - xy (Nx2 array of numbers) – Points positions
- types (Nx1 array of integers) – Type code (1:focus, 2:focus_c, 3:saddle, 4:node_i, 5:node_o)
-
NL_simplifiy
(window_size)[source]¶ Simplify the topological field using Non-local criterions.
Parameters: window_size (number) – Window size used to compute non-local criterions.
-
import_from_CP
(CP_obj, wanted_ind)[source]¶ Import from a CritPoints object, the critical points from the time associated with the given indice.
-
simplify
(dist_min, kind='replacement')[source]¶ Simplify the topological points field.
Parameters: - dist_min (number) – Minimal distance between two points in the simplified field.
- kind (string, optional) – Algorithm used to simplify the field, can be ‘replacement’(default) for iterative replacement with center of mass or ‘only_delete’ to eliminate only the first order associates.
Returns: simpl_TP – Simplified topological points field.
Return type: TopoPoints object
-
class
IMTreatment.vortex_detection.vortex_detection.
VF
(vx, vy, axe_x, axe_y, mask, theta, time)[source]¶ Bases:
object
-
get_cp_cell_position
()[source]¶ Return critical points cell positions and their associated PBI. (PBI : Poincarre_Bendixson indice) Positions are returned in axis unities (axe_x and axe_y) at the center of a cell.
Returns: - pos (2xN array) – position (x, y) of the detected critical points.
- type (1xN array) – type of CP :
ind CP type 0 saddle point 1 unstable focus 2 stable focus 3 unstable node 4 stable node
-
get_cp_position
(thread=1)[source]¶ Return critical points positions and their associated PBI using bilinear interpolation. (PBI : Poincarre_Bendixson indice) Positions are returned in axis unities (axe_x and axe_y).
Parameters: thread (integer or 'all') – Number of thread to use (multiprocessing). Returns: - pos (2xN array) – position (x, y) of the detected critical points.
- pbis (1xN array) – PBI (1 indicate a node, -1 a saddle point)
Note
Using the work of : [1]F. Effenberger and D. Weiskopf, “Finding and classifying critical points of 2D vector fields: a cell-oriented approach using group theory,” Computing and Visualization in Science, vol. 13, no. 8, pp. 377–396, Dec. 2010.
-
-
IMTreatment.vortex_detection.vortex_detection.
get_critical_points
(obj, time=0, unit_time='', window_size=4, kind='pbi', mirroring=None, mirror_interp='linear', smoothing_size=0, verbose=False, thread=1)[source]¶ For a VectorField of a TemporalVectorField object, return the critical points positions and informations on their type.
Parameters: - obj (VectorField or TemporalVectorFields objects) – Base vector field(s)
- time (number, optional) – if ‘obj’ is a VectorField, ‘time’ is the time associated to the field.
- unit_time (string or Unum.units object) – Unit for ‘time’
- window_size (integer, optional) – Size of the interrogation windows for the computation of critical point position (defaut : 4).
- kind (string, optional) – Method used to compute the critical points position. can be : ‘pbi_cell’ for simple (fast) Poincarre-Bendixson sweep. ‘pbi’ for PB sweep and bi-linear interpolation inside cells. ‘pbi_crit’ for PB sweep and use of non-local criterions. ‘gam_vort’ for gamma criterion extremum detection (only detect vortex).
- mirroring (array of numbers) – If specified, use mirroring to get the critical points on the eventual walls. should be an array of ‘[direction (‘x’ or ‘y’), position]*N’.
- mirror_interp (string, optional) – Method used to fill the gap at the wall. ‘value’ : fill with the given value, ‘nearest’ : fill with the nearest value, ‘linear’ (default): fill using linear interpolation (Delaunay triangulation), ‘cubic’ : fill using cubic interpolation (Delaunay triangulation)
- smoothing_size (number, optional) – If specified, a gaussian smoothing of the wanted size is used before detecting the CP.
- verbose (boolean, optional) – If ‘True’, display message on CP detection advancement.
- thread (integer or 'all') – Number of thread to use (multiprocessing). (Only implemented for ‘kind = ‘pbi’‘)
Notes
If the fields have masked values, saddle streamlines ar not computed.
-
IMTreatment.vortex_detection.vortex_detection.
get_vortex_position
(obj, criterion=<function get_residual_vorticity>, criterion_args={}, threshold=0.5, rel=True)[source]¶ Return the position of the vortex (according to the given criterion) on vector field(s).
Parameters: - vectorfield (VectorField or TemporalVectorFields object) –
.
- criterion (function) – Criterion used to highlight vortex position. Should be a function, taking a VectorField object and returning a ScalarField object.
- criterion_args (dict) – Additional arguments to give to the criterion function
- threshold (number) – Threshold value determining the vortex zone.
- rel (Boolean) – If ‘rel’ is ‘True’ (default), ‘threshold’ is relative to the extremum values of the field. If ‘rel’ is ‘False’, ‘threshold’ is treated like an absolut values.
- vectorfield (VectorField or TemporalVectorFields object) –