valib
Vortex Analysis LIBrary
valib_mod.F90
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 
12 module 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_swirling_strength(A, lambda_ci, lambda_cr) &
89  bind(c, name='valib_swirling_strength')
90  use iso_c_binding
91  implicit none
92  real(VA_REAL_F), intent(in) :: A(9)
93  real(VA_REAL_F), intent(out) :: lambda_ci
94  real(VA_REAL_F), intent(out) :: lambda_cr
95  end subroutine valib_swirling_strength
96 
97  subroutine valib_triple_decomposition(A, num_intervals, &
98  residual_vorticity, &
99  residual_strain, &
100  shear) &
101  bind(c, name='valib_triple_decomposition')
102  use iso_c_binding
103  implicit none
104  real(VA_REAL_F), intent(in) :: A(9)
105  integer(c_int), intent(in) :: num_intervals
106  real(VA_REAL_F), intent(out) :: residual_vorticity
107  real(VA_REAL_F), intent(out) :: residual_strain
108  real(VA_REAL_F), intent(out) :: shear
109  end subroutine valib_triple_decomposition
110 
111  subroutine valib_triple_decomposition_4norms(A, num_intervals, &
112  residual_vorticity, &
113  residual_strain, &
114  shear_vorticity, &
115  shear_strain) &
116  bind(c, name='valib_triple_decomposition_4norms')
117  use iso_c_binding
118  implicit none
119  real(VA_REAL_F), intent(in) :: A(9)
120  integer(c_int), intent(in) :: num_intervals
121  real(VA_REAL_F), intent(out) :: residual_vorticity
122  real(VA_REAL_F), intent(out) :: residual_strain
123  real(VA_REAL_F), intent(out) :: shear_vorticity
124  real(VA_REAL_F), intent(out) :: shear_strain
125  end subroutine valib_triple_decomposition_4norms
126  end interface
127 
128 end module valib_mod
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