![]() |
valib
Vortex Analysis LIBrary
|
#include "common_defs.h"
Go to the source code of this file.
Functions | |
| VA_DEVICE_FUN void | valib_contrarotation (VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *S_RAVG) |
| Average contrarotation (a.k.a. \( S_{RAVG} \)) [20]. | |
| VA_DEVICE_FUN void | valib_corotation (VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *omega_RAVG) |
| Average corotation (a.k.a. \(\boldsymbol{\omega}_{RAVG}\)) [13]. | |
| VA_DEVICE_FUN void | valib_delta (VA_REAL *A, VA_REAL *delta) |
| \( \Delta \)-criterion [3] [4] [18] | |
| VA_DEVICE_FUN void | valib_lambda2 (VA_REAL *A, VA_REAL *lambda_2) |
| \( \lambda _2 \)-criterion [9] | |
| VA_DEVICE_FUN void | valib_qcriterion (VA_REAL *A, VA_REAL *qcriterion) |
| \(Q\)-criterion [8] | |
| VA_DEVICE_FUN void | valib_qdcriterion (VA_REAL *A, VA_REAL *qd_criterion) |
| \(Q_D\)-criterion [15] | |
| VA_DEVICE_FUN void | valib_qmcriterion (VA_REAL *A, VA_REAL *qm_criterion) |
| \(Q_M\)-criterion [11] | |
| VA_DEVICE_FUN void | valib_qwcriterion (VA_REAL *A, VA_REAL *qw_criterion) |
| \(Q_W\)-criterion [12] | |
| VA_DEVICE_FUN void | valib_rortex (VA_REAL *A, VA_REAL *rortex) |
| The Rortex vector \( \boldsymbol R \) (a.k.a. the Vortex vector or the Liutex vector) [6] [16]. | |
| VA_DEVICE_FUN void | valib_compactness (VA_REAL *A, VA_REAL *compactness) |
| Compactness criterion. | |
| VA_DEVICE_FUN void | valib_swirling_strength (VA_REAL *A, VA_REAL *lambda_ci, VA_REAL *lambda_cr) |
| Swirling strength criterion (a.k.a. \(\lambda_{ci}\)-criterion) [19] [2]. | |
| VA_DEVICE_FUN void | valib_triple_decomposition (VA_DEVICE_ADDR VA_REAL *A, int *num_intervals, VA_DEVICE_ADDR VA_REAL *residual_vorticity, VA_DEVICE_ADDR VA_REAL *residual_strain, VA_DEVICE_ADDR VA_REAL *shear) |
| Triple Decomposition Method (TDM, a.k.a. residual vorticity and residual strain rate) [14], [10]. | |
| VA_DEVICE_FUN void | valib_residual_vorticity (VA_DEVICE_ADDR VA_REAL *A, int *num_intervals, VA_DEVICE_ADDR VA_REAL *residual_vorticity, VA_DEVICE_ADDR VA_REAL *residual_strain, VA_DEVICE_ADDR VA_REAL *shear) |
| Triple Decomposition Method (TDM, a.k.a. residual vorticity and residual strain rate) [14], [10]. | |
| VA_DEVICE_FUN void | valib_triple_decomposition_4norms (VA_DEVICE_ADDR VA_REAL *A, int *num_intervals, VA_DEVICE_ADDR VA_REAL *residual_vorticity, VA_DEVICE_ADDR VA_REAL *residual_strain, VA_DEVICE_ADDR VA_REAL *shear_vorticity, VA_DEVICE_ADDR VA_REAL *shear_strain) |
| Triple Decomposition Method (TDM, a.k.a. residual vorticity and residual strain rate) [14], [10]. | |
| VA_DEVICE_FUN void | valib_triple_decomposition_with_angles (VA_DEVICE_ADDR VA_REAL *A, int *num_intervals, VA_DEVICE_ADDR VA_REAL *residual_vorticity, VA_DEVICE_ADDR VA_REAL *residual_strain, VA_DEVICE_ADDR VA_REAL *shear, VA_DEVICE_ADDR VA_REAL *aplha, VA_DEVICE_ADDR VA_REAL *beta, VA_DEVICE_ADDR VA_REAL *gamma) |
| Triple Decomposition Method (TDM, a.k.a. residual vorticity and residual strain rate) [14], [10]. | |
Includes for the VALIB library functions.
Definition in file valib.h.