![]() |
valib
Vortex Analysis LIBrary
|
Internal functions used by VALIB. More...
Average corotation and contrarotation helper functions | |
VA_DEVICE_FUN void | valib_resvortstrain2D (VA_REAL *SO, VA_REAL *residual_vorticity, VA_REAL *residual_strain, VA_REAL *principal_axis) |
Numerical integrands on a unit sphere | |
VA_DEVICE_FUN void | valib_integrand_resvort (VA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *avgcorot) |
VA_DEVICE_FUN void | valib_integrand_resstrain (VA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *S_RAVG) |
VA_DEVICE_FUN void | valib_integrand_surface (VA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result) |
Small helper functions | |
VA_DEVICE_FUN VA_REAL | valib_sign (VA_REAL a, VA_REAL b) |
VA_DEVICE_FUN VA_REAL | valib_cubic_root (VA_REAL x) |
VA_DEVICE_FUN VA_REAL | valib_max (VA_REAL a, VA_REAL b) |
VA_DEVICE_FUN VA_REAL | valib_min (VA_REAL a, VA_REAL b) |
VA_DEVICE_FUN void | valib_swap_values (VA_REAL *a, VA_REAL *b) |
VA_DEVICE_FUN void | valib_solve_cubic_equation (VA_REAL a, VA_REAL b, VA_REAL c, VA_REAL d, VA_REAL complex *x1, VA_REAL complex *x2, VA_REAL complex *x3, int *info) |
Operations on vectors of length 3 | |
VA_DEVICE_FUN void | valib_vec_zero3 (VA_REAL *v) |
VA_DEVICE_FUN VA_REAL | valib_scalar_prod3 (VA_REAL *v1, VA_REAL *v2) |
VA_DEVICE_FUN void | valib_cross_prod3 (VA_REAL *v1, VA_REAL *v2, VA_REAL *v3) |
Operations on 3x3 matrices | |
VA_DEVICE_FUN void | valib_mat_zero3 (VA_REAL *A) |
VA_DEVICE_FUN void | valib_mat_copy3 (VA_REAL *A, VA_REAL *B) |
VA_DEVICE_FUN void | valib_norm_frobenius3 (VA_REAL *A, VA_REAL *norm) |
VA_DEVICE_FUN void | valib_trace3 (VA_REAL *A, VA_REAL *tr3) |
VA_DEVICE_FUN void | valib_second_invariant3 (VA_REAL *A, VA_REAL *second_invariant) |
VA_DEVICE_FUN void | valib_determinant3 (VA_REAL *A, VA_REAL *det) |
VA_DEVICE_FUN void | valib_transpose3 (VA_REAL *A) |
VA_DEVICE_FUN void | valib_sym_antisym3 (VA_REAL *A) |
VA_DEVICE_FUN void | valib_deviatorise3 (VA_REAL *A) |
VA_DEVICE_FUN void | valib_matvec3_prod (int trans_a, VA_REAL *A, VA_REAL *x, VA_REAL *y) |
VA_DEVICE_FUN void | valib_rank1_update3 (VA_REAL alpha, VA_REAL *x, VA_REAL *y, VA_REAL *A) |
VA_DEVICE_FUN void | valib_tatt3 (VA_REAL *A, VA_REAL *T, VA_REAL *TATT) |
VA_DEVICE_FUN void | valib_matmat3_prod (int trans_a, int trans_b, VA_REAL *A, VA_REAL *B, VA_REAL *C) |
VA_DEVICE_FUN void | valib_gram_schmidt3 (VA_REAL *Q) |
VA_DEVICE_FUN void | valib_eigenvalues_sym3 (VA_REAL *A, VA_REAL *eigval) |
VA_DEVICE_FUN int | valib_find_pivot3 (VA_REAL *A, int jcol, int start) |
VA_DEVICE_FUN void | valib_find_nullspace3 (VA_REAL *A, VA_REAL *x, int *info) |
VA_DEVICE_FUN void | valib_constructQ3 (VA_REAL sina, VA_REAL cosa, VA_REAL sinb, VA_REAL cosb, VA_REAL sing, VA_REAL cosg, VA_REAL *Q) |
VA_DEVICE_FUN void | valib_constructQfromN3 (VA_REAL *n, VA_REAL *Q) |
Helper functions for the Q-criterion and its variants. | |
VA_DEVICE_FUN void | valib_qcriterion_blending (VA_REAL *SO, VA_REAL blending_ratio, VA_REAL *qcriterion_blended) |
Numerical integrators on a sphere | |
VA_DEVICE_FUN void | valib_integrator_fibonacci (VA_DEVICE_ADDR void integrand(VA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result), VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result) |
VA_DEVICE_FUN void | valib_integrator_bazant (VA_DEVICE_ADDR void integrand(VA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result), VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result) |
VA_DEVICE_FUN void | valib_integrator_rectangle (VA_DEVICE_ADDR void integrand(VA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result), VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result, int numIntervals) |
Triple decomposition helper functions. | |
VA_DEVICE_FUN void | valib_get_brf (VA_DEVICE_ADDR VA_REAL *A, int *num_intervals, VA_REAL *alpha_brf, VA_REAL *beta_brf, VA_REAL *gamma_brf) |
Internal functions used by VALIB.
VA_DEVICE_FUN void valib_resvortstrain2D | ( | VA_REAL * | SO, |
VA_REAL * | residual_vorticity, | ||
VA_REAL * | residual_strain, | ||
VA_REAL * | principal_axis | ||
) |
Function for determining residual vorticity and residual strain in 2D. The 2D case is seen as the leading principal 2x2 block of the 3x3 matrix. The function implements formulas (16)-(27) from [14].
Definition at line 23 of file corotation_common.h.
VA_DEVICE_FUN void valib_integrand_resvort | ( | VA_REAL * | n, |
VA_REAL * | weight, | ||
VA_DEVICE_ADDR VA_REAL * | A, | ||
VA_DEVICE_ADDR VA_REAL * | avgcorot | ||
) |
Integrand to be used with different integration rules for determining averaged corotation in a given integration point with a weight, compute the local contribution and add it to the average corotation vector.
Definition at line 94 of file corotation_common.h.
VA_DEVICE_FUN void valib_integrand_resstrain | ( | VA_REAL * | n, |
VA_REAL * | weight, | ||
VA_DEVICE_ADDR VA_REAL * | A, | ||
VA_DEVICE_ADDR VA_REAL * | S_RAVG | ||
) |
Integrand to be used with different integration rules for determining averaged contrarotation in a given integration point with a weight, compute the local contribution and add it to the average contrarotation.
Definition at line 163 of file corotation_common.h.
VA_DEVICE_FUN void valib_integrand_surface | ( | VA_REAL * | n, |
VA_REAL * | weight, | ||
VA_DEVICE_ADDR VA_REAL * | A, | ||
VA_DEVICE_ADDR VA_REAL * | result | ||
) |
Integrand to be used for debugging purposes. Integrator should return value of the surface of a unit sphere.
Definition at line 269 of file corotation_common.h.
VA_DEVICE_FUN VA_REAL valib_sign | ( | VA_REAL | a, |
VA_REAL | b | ||
) |
VA_DEVICE_FUN VA_REAL valib_cubic_root | ( | VA_REAL | x | ) |
VA_DEVICE_FUN VA_REAL valib_max | ( | VA_REAL | a, |
VA_REAL | b | ||
) |
VA_DEVICE_FUN VA_REAL valib_min | ( | VA_REAL | a, |
VA_REAL | b | ||
) |
VA_DEVICE_FUN void valib_swap_values | ( | VA_REAL * | a, |
VA_REAL * | b | ||
) |
VA_DEVICE_FUN void valib_solve_cubic_equation | ( | VA_REAL | a, |
VA_REAL | b, | ||
VA_REAL | c, | ||
VA_REAL | d, | ||
VA_REAL complex * | x1, | ||
VA_REAL complex * | x2, | ||
VA_REAL complex * | x3, | ||
int * | info | ||
) |
Function for solving cubic equation using Cardano formulas based on K. Rektorys: Prehled uzite matematiky, p. 39 a*x^3 + b*x^2 + c*x + d = 0 If determinant is positive, it cannot solve the problem, info = -1. If determinant is negative, it finds the three roots, info = 0.
VA_DEVICE_FUN void valib_vec_zero3 | ( | VA_REAL * | v | ) |
VA_DEVICE_FUN VA_REAL valib_scalar_prod3 | ( | VA_REAL * | v1, |
VA_REAL * | v2 | ||
) |
VA_DEVICE_FUN void valib_cross_prod3 | ( | VA_REAL * | v1, |
VA_REAL * | v2, | ||
VA_REAL * | v3 | ||
) |
VA_DEVICE_FUN void valib_mat_zero3 | ( | VA_REAL * | A | ) |
VA_DEVICE_FUN void valib_mat_copy3 | ( | VA_REAL * | A, |
VA_REAL * | B | ||
) |
VA_DEVICE_FUN void valib_norm_frobenius3 | ( | VA_REAL * | A, |
VA_REAL * | norm | ||
) |
VA_DEVICE_FUN void valib_trace3 | ( | VA_REAL * | A, |
VA_REAL * | tr3 | ||
) |
VA_DEVICE_FUN void valib_second_invariant3 | ( | VA_REAL * | A, |
VA_REAL * | second_invariant | ||
) |
VA_DEVICE_FUN void valib_determinant3 | ( | VA_REAL * | A, |
VA_REAL * | det | ||
) |
VA_DEVICE_FUN void valib_transpose3 | ( | VA_REAL * | A | ) |
VA_DEVICE_FUN void valib_sym_antisym3 | ( | VA_REAL * | A | ) |
VA_DEVICE_FUN void valib_deviatorise3 | ( | VA_REAL * | A | ) |
VA_DEVICE_FUN void valib_matvec3_prod | ( | int | trans_a, |
VA_REAL * | A, | ||
VA_REAL * | x, | ||
VA_REAL * | y | ||
) |
VA_DEVICE_FUN void valib_rank1_update3 | ( | VA_REAL | alpha, |
VA_REAL * | x, | ||
VA_REAL * | y, | ||
VA_REAL * | A | ||
) |
VA_DEVICE_FUN void valib_tatt3 | ( | VA_REAL * | A, |
VA_REAL * | T, | ||
VA_REAL * | TATT | ||
) |
VA_DEVICE_FUN void valib_matmat3_prod | ( | int | trans_a, |
int | trans_b, | ||
VA_REAL * | A, | ||
VA_REAL * | B, | ||
VA_REAL * | C | ||
) |
VA_DEVICE_FUN void valib_gram_schmidt3 | ( | VA_REAL * | Q | ) |
VA_DEVICE_FUN void valib_eigenvalues_sym3 | ( | VA_REAL * | A, |
VA_REAL * | eigval | ||
) |
VA_DEVICE_FUN int valib_find_pivot3 | ( | VA_REAL * | A, |
int | jcol, | ||
int | start | ||
) |
VA_DEVICE_FUN void valib_find_nullspace3 | ( | VA_REAL * | A, |
VA_REAL * | x, | ||
int * | info | ||
) |
VA_DEVICE_FUN void valib_constructQ3 | ( | VA_REAL | sina, |
VA_REAL | cosa, | ||
VA_REAL | sinb, | ||
VA_REAL | cosb, | ||
VA_REAL | sing, | ||
VA_REAL | cosg, | ||
VA_REAL * | Q | ||
) |
VA_DEVICE_FUN void valib_constructQfromN3 | ( | VA_REAL * | n, |
VA_REAL * | Q | ||
) |
|
inline |
Function common for \(Q\), \(Q_D\), \(Q_W\), and \(Q_M\) criterions with prescribed blending factor \(\beta\) of the norms \( \| \Omega \|_F^2 \) and \( \| S_X \|_F^2\). Evaluates
\[ Q = \frac{1}{2} ( \| \Omega \|_F^2 - \beta \| S_X \|_F^2 ) > 0, \]
where \( \Omega = \frac{1}{2} ( \nabla u - ( \nabla u )^T ) \) is the antisymmetric part and \( S = \frac{1}{2} ( \nabla u + ( \nabla u )^T ) \) is the symmetric part of the velocity gradient \( \nabla u \).
For the \(Q\)-criterion [8] \( S_X = S \) and \(\beta = 1\).
For the \(Q_D\)-criterion [15] \( S_X = S_D \) and \(\beta = 1\).
For the \(Q_M\)-criterion [11] \( S_X = S_D \) and \(\beta = 1.5\).
For the \(Q_W\)-criterion [12] \( S_X = S_D \) and \(\beta = 1.2\).
Here, \( S_D \) is the deviatoric part of \( S \), which is the same as the symmetric part of the deviatoric part of the velocity gradient,
\[ ( \nabla u )_D = \nabla u - \frac{1}{3} \mbox{trace}(\nabla u)\ I. \]
Definition at line 50 of file qcriterion_common.h.
VA_DEVICE_FUN void valib_integrator_fibonacci | ( | VA_DEVICE_ADDR void | integrandVA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result, |
VA_DEVICE_ADDR VA_REAL * | A, | ||
VA_DEVICE_ADDR VA_REAL * | result | ||
) |
Integrator on a unit sphere using Fibonacci numbers, [7]. Programmed following Appendix A from [5].
Definition at line 20 of file spherical_integrators.h.
VA_DEVICE_FUN void valib_integrator_bazant | ( | VA_DEVICE_ADDR void | integrandVA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result, |
VA_DEVICE_ADDR VA_REAL * | A, | ||
VA_DEVICE_ADDR VA_REAL * | result | ||
) |
Integrator on a unit sphere following [1].
Definition at line 74 of file spherical_integrators.h.
VA_DEVICE_FUN void valib_integrator_rectangle | ( | VA_DEVICE_ADDR void | integrandVA_REAL *n, VA_REAL *weight, VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *result, |
VA_DEVICE_ADDR VA_REAL * | A, | ||
VA_DEVICE_ADDR VA_REAL * | result, | ||
int | numIntervals | ||
) |
Integrator on a unit sphere using the simplest rectangular rule. Corresponds to the approach described in [13].
Definition at line 116 of file spherical_integrators.h.
VA_DEVICE_FUN void valib_get_brf | ( | VA_DEVICE_ADDR VA_REAL * | A, |
int * | num_intervals, | ||
VA_REAL * | alpha_brf, | ||
VA_REAL * | beta_brf, | ||
VA_REAL * | gamma_brf | ||
) |
Function for determining the Basic Reference Frame (BRF) introduced in [14]. This is the coordinate frame in which the shear is maximized and in which it is later eliminated. The BRF is found as an optimization problem, searching the maxima in formula (10) of [14], i.e.
\[ \max_{\alpha \in [0,\pi], \\ \beta \in [0,\pi], \\ \gamma \in [0,\pi/2]} \left( |S_{12}\Omega_{12}| + |S_{23}\Omega_{23}| + |S_{31}\Omega_{31}| \right), \]
where \(\alpha\), \(\beta\), and \(\gamma\) are the three angles of rotation of the coordinate frame. The optimal value is determined by an exhaustive search over all possible values of the angles. The resolution is determined as \( h = \pi / n_i \), with \(n_i\) the number of intervals (the num_intervals variable). For example, num_intervals = 180 corresponds to 1 degree resolution. This optimization problem can become costly for large number of intervals.
Definition at line 38 of file triple_decomposition.h.