valib
Vortex Analysis LIBrary
Loading...
Searching...
No Matches
valib_mod.F90
Go to the documentation of this file.
1
3
4!*******************************************************************************
5! Edit to switch to single-precision version (double is the default):
6! double precision version
7!#define VA_REAL_F c_double
8! single precision version
9#define VA_REAL_F c_float
10!*******************************************************************************
11
12module valib_mod
13 use iso_c_binding
14 implicit none
15 interface
16 subroutine valib_contrarotation(A, S_RAVG) &
17 bind(c, name='valib_contrarotation')
18 use iso_c_binding
19 implicit none
20 real(VA_REAL_F), intent(in) :: A(9)
21 real(VA_REAL_F), intent(out) :: S_RAVG(9)
22 end subroutine valib_contrarotation
23
24 subroutine valib_corotation(A, omega_RAVG) &
25 bind(c, name='valib_corotation')
26 use iso_c_binding
27 implicit none
28 real(VA_REAL_F), intent(in) :: A(9)
29 real(VA_REAL_F), intent(out) :: omega_RAVG(3)
30 end subroutine valib_corotation
31
32 subroutine valib_delta(A, delta) &
33 bind(c, name='valib_delta')
34 use iso_c_binding
35 implicit none
36 real(VA_REAL_F), intent(in) :: A(9)
37 real(VA_REAL_F), intent(out) :: delta
38 end subroutine valib_delta
39
40 subroutine valib_lambda2(A, lambda_2) &
41 bind(c, name='valib_lambda2')
42 use iso_c_binding
43 implicit none
44 real(VA_REAL_F), intent(in) :: A(9)
45 real(VA_REAL_F), intent(out) :: lambda_2
46 end subroutine valib_lambda2
47
48 subroutine valib_qcriterion(A, qcriterion) &
49 bind(c, name='valib_qcriterion')
50 use iso_c_binding
51 implicit none
52 real(VA_REAL_F), intent(in) :: A(9)
53 real(VA_REAL_F), intent(out) :: qcriterion
54 end subroutine valib_qcriterion
55
56 subroutine valib_qdcriterion(A, qdcriterion) &
57 bind(c, name='valib_qdcriterion')
58 use iso_c_binding
59 implicit none
60 real(VA_REAL_F), intent(in) :: A(9)
61 real(VA_REAL_F), intent(out) :: qdcriterion
62 end subroutine valib_qdcriterion
63
64 subroutine valib_qmcriterion(A, qmcriterion) &
65 bind(c, name='valib_qmcriterion')
66 use iso_c_binding
67 implicit none
68 real(VA_REAL_F), intent(in) :: A(9)
69 real(VA_REAL_F), intent(out) :: qmcriterion
70 end subroutine valib_qmcriterion
71
72 subroutine valib_qwcriterion(A, qwcriterion) &
73 bind(c, name='valib_qwcriterion')
74 use iso_c_binding
75 implicit none
76 real(VA_REAL_F), intent(in) :: A(9)
77 real(VA_REAL_F), intent(out) :: qwcriterion
78 end subroutine valib_qwcriterion
79
80 subroutine valib_rortex(A, rortex) &
81 bind(c, name='valib_rortex')
82 use iso_c_binding
83 implicit none
84 real(VA_REAL_F), intent(in) :: A(9)
85 real(VA_REAL_F), intent(out) :: rortex(3)
86 end subroutine valib_rortex
87
88 subroutine valib_compactness(A, compactness) &
89 bind(c, name='valib_compactness')
90 use iso_c_binding
91 implicit none
92 real(VA_REAL_F), intent(in) :: A(9)
93 real(VA_REAL_F), intent(out) :: compactness
94 end subroutine valib_compactness
95
96 subroutine valib_swirling_strength(A, lambda_ci, lambda_cr) &
97 bind(c, name='valib_swirling_strength')
98 use iso_c_binding
99 implicit none
100 real(VA_REAL_F), intent(in) :: A(9)
101 real(VA_REAL_F), intent(out) :: lambda_ci
102 real(VA_REAL_F), intent(out) :: lambda_cr
103 end subroutine valib_swirling_strength
104
105 subroutine valib_triple_decomposition(A, num_intervals, &
106 residual_vorticity, &
107 residual_strain, &
108 shear) &
109 bind(c, name='valib_triple_decomposition')
110 use iso_c_binding
111 implicit none
112 real(VA_REAL_F), intent(in) :: A(9)
113 integer(c_int), intent(in) :: num_intervals
114 real(VA_REAL_F), intent(out) :: residual_vorticity
115 real(VA_REAL_F), intent(out) :: residual_strain
116 real(VA_REAL_F), intent(out) :: shear
117 end subroutine valib_triple_decomposition
118
119 subroutine valib_triple_decomposition_4norms(A, num_intervals, &
120 residual_vorticity, &
121 residual_strain, &
122 shear_vorticity, &
123 shear_strain) &
124 bind(c, name='valib_triple_decomposition_4norms')
125 use iso_c_binding
126 implicit none
127 real(VA_REAL_F), intent(in) :: A(9)
128 integer(c_int), intent(in) :: num_intervals
129 real(VA_REAL_F), intent(out) :: residual_vorticity
130 real(VA_REAL_F), intent(out) :: residual_strain
131 real(VA_REAL_F), intent(out) :: shear_vorticity
132 real(VA_REAL_F), intent(out) :: shear_strain
134
135 subroutine valib_triple_decomposition_with_angles(A, num_intervals, &
136 residual_vorticity, &
137 residual_strain, &
138 shear, &
139 alpha, &
140 beta, &
141 gamma) &
142 bind(c, name='valib_triple_decomposition_with_angles')
143 use iso_c_binding
144 implicit none
145 real(VA_REAL_F), intent(in) :: A(9)
146 integer(c_int), intent(in) :: num_intervals
147 real(VA_REAL_F), intent(out) :: residual_vorticity
148 real(VA_REAL_F), intent(out) :: residual_strain
149 real(VA_REAL_F), intent(out) :: shear
150 real(VA_REAL_F), intent(out) :: alpha
151 real(VA_REAL_F), intent(out) :: beta
152 real(VA_REAL_F), intent(out) :: gamma
154 end interface
155
156end module valib_mod
VA_DEVICE_FUN void valib_compactness(VA_REAL *A, VA_REAL *compactness)
Compactness criterion.
Definition compactness.h:38
VA_DEVICE_FUN void valib_contrarotation(VA_DEVICE_ADDR VA_REAL *A, VA_DEVICE_ADDR VA_REAL *S_RAVG)
Average contrarotation (a.k.a. ) .
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) ,...
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
VA_DEVICE_FUN void valib_qwcriterion(VA_REAL *A, VA_REAL *qw_criterion)
-criterion
Definition qwcriterion.h:48
VA_DEVICE_FUN void valib_delta(VA_REAL *A, VA_REAL *delta)
-criterion
Definition delta.h:49
VA_DEVICE_FUN void valib_lambda2(VA_REAL *A, VA_REAL *lambda_2)
-criterion
Definition lambda2.h:59
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) ,...
VA_DEVICE_FUN void valib_qcriterion(VA_REAL *A, VA_REAL *qcriterion)
-criterion
Definition qcriterion.h:46
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) .
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 *alpha, 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) ,...
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
VA_DEVICE_FUN void valib_qdcriterion(VA_REAL *A, VA_REAL *qd_criterion)
-criterion
Definition qdcriterion.h:45
VA_DEVICE_FUN void valib_qmcriterion(VA_REAL *A, VA_REAL *qm_criterion)
-criterion
Definition qmcriterion.h:54