numbacs.diagnostics
Module Contents
- numbacs.diagnostics.ftle_grid_2D(flowmap, T, dx, dy)[source]
Compute 2D FTLE field from flowmap which is solution of ode over an initial grid defined by x and y for integration time T.
- Parameters:
flowmap (np.ndarray, shape = (nx,ny,2)) – array containing final positions of trajectories of initial grid from t0 to t0+T.
T (float) – integration time.
dx (float) – grid spacing in x-direction.
dy (float) – grid spacing in y-direction.
- Returns:
ftle – array containing ftle values.
- Return type:
np.ndarray, shape = (nx,ny)
- numbacs.diagnostics.C_tensor_2D(flowmap_aux, dx, dy, h=1e-05)[source]
Compute eigenvalues and eigenvectors of Cauchy Green tensor in 2D from flowmap_aux which is solution of ode over an auxilary grid defined by x,y +-h for integration time T.
- Parameters:
flowmap_aux (np.ndarray, shape = (nx,ny,n_aux,2)) – array containing final positions of initial grid with aux grid spacing h from t0 to t0+T.
dx (float) – grid spacing in x-direction.
dy (float) – grid spacing in y-direction.
h (float, optional) – aux grid spacing. The default is 1e-5.
- Returns:
C – array containing C11, C12, C22 components of Cauchy Green tensor.
- Return type:
np.ndarray, shape = (nx,ny,3)
- numbacs.diagnostics.C_eig_aux_2D(flowmap_aux, dx, dy, h=1e-05, eig_main=True)[source]
Compute eigenvalues and eigenvectors of Cauchy Green tensor in 2D from flowmap_aux which is solution of ode over an auxilary grid defined by x,y +-h for integration time T.
- Parameters:
flowmap_aux (np.ndarray, shape = (nx,ny,n_aux,2)) – array containing final positions of initial grid with aux grid spacing h from t0 to t0+T.
dx (float) – grid spacing in x-direction.
dy (float) – grid spacing in y-direction.
h (float, optional) – aux grid spacing. The default is 1e-5.
eig_main (boolean, optional) – flag to determine if eigevalues are computed from main grid. The default is True.
- Returns:
eigvals (np.ndarray, shape = (nx,ny,2)) – array containing eigenvalues values.
eigvecs (np.ndarray, shape = (nx,ny,2,2)) – array containing eigenvectors.
- numbacs.diagnostics.C_eig_2D(flowmap, dx, dy)[source]
Compute eigenvalues and eigenvectors of Cauchy Green tensor in 2D from flowmap which is solution of ode over a grid defined by x,y for integration time T.
- Parameters:
flowmap (np.ndarray, shape = (nx,ny,2)) – array containing final positions of initial grid from t0 to t0+T.
dx (float) – grid spacing in x-direction.
dy (float) – grid spacing in y-direction.
- Returns:
eigvals (np.ndarray, shape = (nx,ny,2)) – array containing eigenvalues.
eigvecs (np.ndarray, shape = (nx,ny,2,2)) – array containing eigenvectors.
- numbacs.diagnostics.ftle_from_eig(eigval_max, T)[source]
Compute FTLE field from eigval_max where eigval_max is eigenvalue of Cauchy-Green tensor computed using integration time T.
- Parameters:
eigval_max (np.ndarray, shape = (nx,ny)) – maximum eigenvalue of Cauchy-Green tensor computed using integration time T.
T (float) – integration time.
- Returns:
ftle – array containing ftle values.
- Return type:
np.ndarray, shape = (nx,ny)
- numbacs.diagnostics.lavd_grid_2D(flowmap_n, tspan, T, vort_interp, xrav, yrav, period_x=0.0, period_y=0.0)[source]
Compute LAVD from flowmap_n where flowmap_n contains trajectories computed over gridpoints defined by xrav,yrav for an integration time T and trajectories are returned at times given by tspan. vort_interp is an interpolant function of vorticity over (at least) that time window.
- Parameters:
flowmap_n (np.ndarray, shape = (nx,ny,n,2)) – array containing trajectories of initial grid from t0 to t0+T.
tspan (np.ndarray, shape = (n,)) – array containing n times corresponding to axis=2 of flowmap_n.
T (float) – integration time.
vort_interp (jit-callable) – interpolant function of vorticity which must be (at least) defined over all values of flowmap_n and times from tspan.
xrav (np.ndarray, shape = (nx*ny,))) – array containing raveled or flattened meshgrid X, can be obtained using X.ravel() or X.flatten().
yrav (np.ndarray, shape = (nx*ny,))) – array containing raveled or flattened meshgrid Y, can be obtained using Y.ravel() or Y.flatten().
period_x (float) – value for period in x-direction, if not periodic, set equal to 0.0. The default is 0.0.
period_y (float) – value for period in y-direction, if not periodic, set equal to 0.0. The default is 0.0.
- Returns:
lavd – array containing lavd values.
- Return type:
np.ndarray, shape = (nx,ny)
- numbacs.diagnostics.ftle_grid_ND(flowmap, IC, T, dX)[source]
Compute ND FTLE field from flowmap which is solution of ode over an initial grid defined by IC for integration time T.
- Parameters:
flowmap (np.ndarray, shape = (nx_1*nx_2*...*nx_ndims,ndims)) – array containing final positions of initial grid from t0 to t0+T.
IC (np.ndarray, shape = (nx_1,nx_2,...,nx_ndims,ndims)) – initial condition array.
T (float) – integration time.
dX (np.ndarray, shape = (ndims,)) – array containing spacing in each x_i direction for i = 1,…,ndims.
- Returns:
ftle – array containing ftle values.
- Return type:
np.ndarray, shape = (nx_1*nx_2*…*nx_ndims,)
- numbacs.diagnostics.ile_2D_func(vel, x, y, t0=None, h=0.001)[source]
Compute the iLE field from the flow defined by vel which is a jit-callable function, step size of h is used in finite differencing.
- Parameters:
vel (jit-callable) – function returing interpolated of function value of velocity in x,y directions.
x (np.ndarray, shape = (nx,)) – array containing x-values.
y (np.ndarray, shape = (ny,)) – array containing y-values.
t0 (float or None, optional) – time value at which to evaluate if v1,v2 interpolants depend on time, if they do not depend on time, set to None. The default is None.
h (float, optional) – step size to be used in finite differencing. The default is 1e-3.
- Returns:
ile – array containing ile values.
- Return type:
np.ndarray, shape = (nx,ny)
- numbacs.diagnostics.S_eig_2D_func(vel, x, y, t0=None, h=0.001)[source]
Compute eigenvalues and eigenvectors of Eulerian rate of strain tensor in 2D from vel which is a jit-callable function, step size of h is used for finite differencing.
- Parameters:
vel (jit-callable) – function returing interpolated of function value of velocity in x,y directions.
x (np.ndarray, shape = (nx,)) – array containing x-values.
y (np.ndarray, shape = (ny,)) – array containing y-values.
t0 (float or None, optional) – time value at which to evaluate if v1,v2 interpolants depend on time, if they do not depend on time, set to None. The default is None.
h (float, optional) – step size to be used in finite differencing. The default is 1e-2.
- Returns:
eigvals (np.ndarray, shape = (nx,ny)) – array containg eigenvalues of S.
eigvecs (np.ndarray, shape = (nx,ny)) – array containing eigenvectors of S.
- numbacs.diagnostics.S_2D_func(vel, x, y, t0=None, h=0.001)[source]
Compute Eulerian rate of strain tensor in 2D from vel which is a jit-callable functions, step size of h is used for finite differencing.
- Parameters:
v1 (jit-callable) – function returing interpolated value of velocity in x-direction.
v2 (jit-callable) – function returing interpolated value of velocity in y-direction.
x (np.ndarray, shape = (nx,)) – array containing x-values.
y (np.ndarray, shape = (ny,)) – array containing y-values.
t0 (float or None, optional) – time value at which to evaluate if v1,v2 interpolants depend on time, if they do not depend on time, set to None. The default is None.
h (float, optional) – step size to be used in finite differencing. The default is 1e-3.
- Returns:
S – S11,S12, and S22 components of S tensor.
- Return type:
np.ndarray, shape = (nx,ny)
- numbacs.diagnostics.ile_2D_data(u, v, dx, dy)[source]
Compute the iLE field from the flow defined by u,v over an intial grid defined by x,y.
- Parameters:
u (np.ndarray, shape = (nx,ny)) – array containing velocity values in x-direction.
v (np.ndarray, shape = (nx,ny)) – array containing velocity values in y-direction.
dx (float) – grid spacing in x-direction.
dy (float) – grid spacing in y-direction.
- Returns:
ile – array containing ile values.
- Return type:
np.ndarray, shape = (nx,ny)
- numbacs.diagnostics.S_eig_2D_data(u, v, dx, dy)[source]
Compute eigenvalues and eigenvectors of Eulerian rate of strain tensor in 2D from u,v which describe the x,y velocity.
- Parameters:
u (np.ndarray, shape = (nx,ny)) – array containing velocity values in x-direction.
v (np.ndarray, shape = (nx,ny)) – array containing velocity values in y-direction.
dx (float) – grid spacing in x-direction.
dy (float) – grid spacing in y-direction.
- Returns:
eigvals (np.ndarray, shape = (nx,ny,2)) – array containg eigenvalues of S.
eigvecs (np.ndarray, shape = (nx,ny,2,2)) – array containing eigenvectors of S.
- numbacs.diagnostics.ivd_grid_2D(vort, vort_avg)[source]
Compute IVD at an instant from vort and vort_avg.
- Parameters:
vort (np.ndarray, shape = (nx,ny)) – array containg vorticity values.
vort_avg (float) – value of average spatial vorticity at specific time.
- Returns:
ivd – array containing values of ivd.
- Return type:
np.ndarray, shape = (nx,ny)