valib
Vortex Analysis LIBrary
valib.h
1 
7 #ifndef VA_VALIB_H
8 #define VA_VALIB_H
9 
10 #include "common_defs.h"
11 
12 VA_DEVICE_FUN void valib_contrarotation(VA_DEVICE_ADDR VA_REAL *A,
13  VA_DEVICE_ADDR VA_REAL *S_RAVG);
14 
15 VA_DEVICE_FUN void valib_corotation(VA_DEVICE_ADDR VA_REAL *A,
16  VA_DEVICE_ADDR VA_REAL *omega_RAVG);
17 
18 VA_DEVICE_FUN void valib_delta(VA_REAL *A, VA_REAL *delta);
19 
20 VA_DEVICE_FUN void valib_lambda2(VA_REAL *A, VA_REAL *lambda_2);
21 
22 VA_DEVICE_FUN void valib_qcriterion(VA_REAL *A, VA_REAL *qcriterion);
23 
24 VA_DEVICE_FUN void valib_qdcriterion(VA_REAL *A, VA_REAL *qd_criterion);
25 
26 VA_DEVICE_FUN void valib_qmcriterion(VA_REAL *A, VA_REAL *qm_criterion);
27 
28 VA_DEVICE_FUN void valib_qwcriterion(VA_REAL *A, VA_REAL *qw_criterion);
29 
30 VA_DEVICE_FUN void valib_rortex(VA_REAL *A, VA_REAL *rortex);
31 
32 VA_DEVICE_FUN void valib_swirling_strength(VA_REAL *A,
33  VA_REAL *lambda_ci,
34  VA_REAL *lambda_cr);
35 
36 VA_DEVICE_FUN void valib_triple_decomposition(
37  VA_DEVICE_ADDR VA_REAL *A,
38  int *num_intervals,
39  VA_DEVICE_ADDR VA_REAL *residual_vorticity,
40  VA_DEVICE_ADDR VA_REAL *residual_strain,
41  VA_DEVICE_ADDR VA_REAL *shear);
42 
43 VA_DEVICE_FUN void valib_triple_decomposition_4norms(
44  VA_DEVICE_ADDR VA_REAL *A, int *num_intervals,
45  VA_DEVICE_ADDR VA_REAL *residual_vorticity,
46  VA_DEVICE_ADDR VA_REAL *residual_strain,
47  VA_DEVICE_ADDR VA_REAL *shear_vorticity,
48  VA_DEVICE_ADDR VA_REAL *shear_strain);
49 
50 #endif // VA_VALIB_H
valib_corotation
VA_DEVICE_FUN void valib_corotation(VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *omega_RAVG)
Average corotation (a.k.a. ) .
Definition: corotation.h:89
valib_triple_decomposition
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) ,...
Definition: triple_decomposition.h:193
valib_qwcriterion
VA_DEVICE_FUN void valib_qwcriterion(VA_REAL *A, VA_REAL *qw_criterion)
-criterion
Definition: qwcriterion.h:48
valib_contrarotation
VA_DEVICE_FUN void valib_contrarotation(VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *S_RAVG)
Average contrarotation (a.k.a. ) .
Definition: contrarotation.h:55
valib_qmcriterion
VA_DEVICE_FUN void valib_qmcriterion(VA_REAL *A, VA_REAL *qm_criterion)
-criterion
Definition: qmcriterion.h:54
valib_qcriterion
VA_DEVICE_FUN void valib_qcriterion(VA_REAL *A, VA_REAL *qcriterion)
-criterion
Definition: qcriterion.h:46
valib_triple_decomposition_4norms
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) ,...
Definition: triple_decomposition.h:340
valib_swirling_strength
VA_DEVICE_FUN void valib_swirling_strength(VA_REAL *A, VA_REAL *lambda_ci, VA_REAL *lambda_cr)
Swirling strength criterion (a.k.a. -criterion) .
Definition: swirling_strength.h:56
valib_qdcriterion
VA_DEVICE_FUN void valib_qdcriterion(VA_REAL *A, VA_REAL *qd_criterion)
-criterion
Definition: qdcriterion.h:45
valib_rortex
VA_DEVICE_FUN void valib_rortex(VA_REAL *A, VA_REAL *rortex)
The Rortex vector (a.k.a. the Vortex vector or the Liutex vector) .
Definition: rortex.h:53
valib_delta
VA_DEVICE_FUN void valib_delta(VA_REAL *A, VA_REAL *delta)
-criterion
Definition: delta.h:49
valib_lambda2
VA_DEVICE_FUN void valib_lambda2(VA_REAL *A, VA_REAL *lambda_2)
-criterion
Definition: lambda2.h:59