viscid.calculator.minvar_tools module

Minimum Variance Analysis and boundary normal crd tools

viscid.calculator.minvar_tools.minvar(B, p1, p2, n=40)[source]

Find minimum variance eigenvectors of B between p1 and p2

Do minimum variance analysis (MVA) of vector field B between two points, i.e., p1 and p2 are like two spacecraft locations on either side of the boundary for the MVA.

Parameters:
  • B (viscid.field.VectorField) – Vector field for MVA
  • p1 (sequence) – xyz point on one side of the boundary
  • p2 (sequence) – xyz point on the other side of the boundary
  • n (int) – number of points to sample B for the MVA
Returns:

(evals, evecs) All are ordered toward increasing eigenvalue magnitude. evecs is a 2d array where the vectors are columns (2nd axis), i.e., the min variance eigenvector is evecs[0, :]

viscid.calculator.minvar_tools.minvar_series(arr, warn_finite=True)[source]

Find minimum variance eigenvectors of vector series

Parameters:arr (ndarray) – vector time series with shape (nsamples, 3)
Returns:(evals, evecs) All are ordered toward increasing eigenvalue magnitude. evecs is a 2d array where the vectors are columns (2nd axis), i.e., the min variance eigenvector is evecs[0, :]
viscid.calculator.minvar_tools.minvar_around(B, p0, l=1.0, path_dir=(1, 0, 0), n=40)[source]

Find minimum variance eigenvectors of B around point p0

Do minimum variance analysis (MVA) of vector field B a distance l/2 around point p0 in the direction of path_dir, i.e., path_dir is like the spacecraft travel direction.

Parameters:
  • B (viscid.field.VectorField) – Vector field for MVA
  • p0 (sequence) – xyz center point, within -l/2 and +l/2 of the boundary in the path_dir directon
  • l (float) – distance on either side of center point
  • path_dir (tuple) – spacecraft travel direction
  • n (int) – number of points to sample B for the MVA
Returns:

(evals, evecs) All are ordered toward increasing eigenvalue magnitude. evecs is a 2d array where the vectors are columns (2nd axis), i.e., the min variance eigenvector is evecs[0, :]

See also

minvar()

viscid.calculator.minvar_tools.find_minvar_lmn(B, p1, p2, n=40, l_basis=(0, 0, 1))[source]

Find rotation matrix for going to boundary normal (lmn) crds

Parameters:
  • B (viscid.field.VectorField) – Vector field for MVA
  • p1 (sequence) – xyz point on one side of the boundary
  • p2 (sequence) – xyz point on the other side of the boundary
  • n (int) – number of points to sample B for the MVA
  • l_basis (sequence) – vector used to define the l-direction, i.e., the l direction will be l_basis projected perpendicular to the boundary normal directon. If None, then l_basis will be set to evec_max.
Returns:

3x3 array of l,m,n basis column (2nd axis) vectors

Return type:

ndarray

viscid.calculator.minvar_tools.find_minvar_lmn_around(B, p0, l=1.0, path_dir=(1, 0, 0), n=40, l_basis=(0, 0, 1))[source]

Find rotation matrix for going to boundary normal (lmn) crds

Finds boundary normal from minimum variance analysis (MVA) of vector field B a distance l/2 around point p0 in the direction of path_dir, i.e., path_dir is like the spacecraft travel direction.

Parameters:
  • B (viscid.field.VectorField) – Vector field for MVA
  • p0 (sequence) – xyz center point, within -l/2 and +l/2 of the boundary in the path_dir directon
  • l (float) – distance on either side of center point
  • path_dir (tuple) – spacecraft travel direction
  • n (int) – number of points to sample B for the MVA
  • l_basis (sequence) – vector used to define the l-direction, i.e., the l direction will be l_basis projected perpendicular to the boundary normal directon. If None, then l_basis will be set to evec_max.
Returns:

3x3 array of l,m,n basis column (2nd axis) vectors

Return type:

ndarray

See also

find_nlm()