Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
vreman_cpu.f90
Go to the documentation of this file.
1! Copyright (c) 2023-2026, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use num_types, only : rp
36 use field_list, only : field_list_t
37 use math, only : cadd, neko_eps, col2, vlsc2
39 use registry, only : neko_registry
40 use field, only : field_t
41 use operators, only : dudxyz, grad
42 use coefs, only : coef_t
43 use gs_ops, only : gs_op_add
44 use utils, only : neko_error
45 implicit none
46 private
47
48 public :: vreman_compute_cpu
49
50contains
51
65 subroutine vreman_compute_cpu(if_ext, t, tstep, coef, nut, delta, c, &
66 if_corr, scalar_name, ri_c, ref_temp, g)
67 logical, intent(in) :: if_ext, if_corr
68 real(kind=rp), intent(in) :: t
69 integer, intent(in) :: tstep
70 type(coef_t), intent(in) :: coef
71 type(field_t), intent(inout) :: nut
72 type(field_t), intent(in) :: delta
73 character(len=*), intent(in) :: scalar_name
74 real(kind=rp), intent(in) :: c, ri_c, ref_temp
75 real(kind=rp), intent(in) :: g(3)
76 ! This is the alpha tensor in the paper
77 type(field_t), pointer :: a11, a12, a13, a21, a22, a23, a31, a32, a33
78 type(field_t), pointer :: u, v, w
79 type(field_t), pointer :: temperature, dtdx, dtdy, dtdz
80
81 real(kind=rp) :: beta11
82 real(kind=rp) :: beta12
83 real(kind=rp) :: beta13
84 real(kind=rp) :: beta22
85 real(kind=rp) :: beta23
86 real(kind=rp) :: beta33
87 real(kind=rp) :: b_beta
88 real(kind=rp) :: aijaij
89 integer :: temp_indices(9)
90 integer :: temp_indices_buoy(3)
91 integer :: e, i, j
92 real(kind=rp) :: gmag, ri, correction, buoyancy, shear_sq
93 real(kind=rp) :: n(3), du_n(3), sh(3)
94 real(kind=rp) :: du_parallel
95
96 if (if_ext) then
97 u => neko_registry%get_field_by_name("u_e")
98 v => neko_registry%get_field_by_name("v_e")
99 w => neko_registry%get_field_by_name("w_e")
100 else
101 u => neko_registry%get_field_by_name("u")
102 v => neko_registry%get_field_by_name("v")
103 w => neko_registry%get_field_by_name("w")
104 end if
105
106 call neko_scratch_registry%request_field(a11, temp_indices(1), .false.)
107 call neko_scratch_registry%request_field(a12, temp_indices(2), .false.)
108 call neko_scratch_registry%request_field(a13, temp_indices(3), .false.)
109 call neko_scratch_registry%request_field(a21, temp_indices(4), .false.)
110 call neko_scratch_registry%request_field(a22, temp_indices(5), .false.)
111 call neko_scratch_registry%request_field(a23, temp_indices(6), .false.)
112 call neko_scratch_registry%request_field(a31, temp_indices(7), .false.)
113 call neko_scratch_registry%request_field(a32, temp_indices(8), .false.)
114 call neko_scratch_registry%request_field(a33, temp_indices(9), .false.)
115
116
117 ! Compute the derivatives of the velocity (the alpha tensor)
118 call dudxyz (a11%x, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
119 call dudxyz (a12%x, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
120 call dudxyz (a13%x, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
121
122 call dudxyz (a21%x, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
123 call dudxyz (a22%x, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
124 call dudxyz (a23%x, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
125
126 call dudxyz (a31%x, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
127 call dudxyz (a32%x, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
128 call dudxyz (a33%x, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
129
130 call coef%gs_h%op(a11, gs_op_add)
131 call coef%gs_h%op(a12, gs_op_add)
132 call coef%gs_h%op(a13, gs_op_add)
133 call coef%gs_h%op(a21, gs_op_add)
134 call coef%gs_h%op(a22, gs_op_add)
135 call coef%gs_h%op(a23, gs_op_add)
136 call coef%gs_h%op(a31, gs_op_add)
137 call coef%gs_h%op(a32, gs_op_add)
138 call coef%gs_h%op(a33, gs_op_add)
139
140 !$omp parallel do private(e, i, beta11, beta12, beta13, beta22, beta23, &
141 !$omp& beta33, b_beta, aijaij)
142 do e = 1, coef%msh%nelv
143 !OCL NORECURRENCE, NOVREC, NOALIAS
144 !DIR$ CONCURRENT
145 !GCC$ ivdep
146 do i = 1, coef%Xh%lxyz
147 ! beta_ij = alpha_mi alpha_mj
148 beta11 = a11%x(i,1,1,e)**2 + a21%x(i,1,1,e)**2 + a31%x(i,1,1,e)**2
149 beta22 = a12%x(i,1,1,e)**2 + a22%x(i,1,1,e)**2 + a32%x(i,1,1,e)**2
150 beta33 = a13%x(i,1,1,e)**2 + a23%x(i,1,1,e)**2 + a33%x(i,1,1,e)**2
151 beta12 = a11%x(i,1,1,e)*a12%x(i,1,1,e) + &
152 a21%x(i,1,1,e)*a22%x(i,1,1,e) + &
153 a31%x(i,1,1,e)*a32%x(i,1,1,e)
154 beta13 = a11%x(i,1,1,e)*a13%x(i,1,1,e) + &
155 a21%x(i,1,1,e)*a23%x(i,1,1,e) + &
156 a31%x(i,1,1,e)*a33%x(i,1,1,e)
157 beta23 = a12%x(i,1,1,e)*a13%x(i,1,1,e) + &
158 a22%x(i,1,1,e)*a23%x(i,1,1,e) + &
159 a32%x(i,1,1,e)*a33%x(i,1,1,e)
160
161 b_beta = beta11*beta22 - beta12*beta12 + beta11*beta33 &
162 - beta13*beta13 + beta22*beta33 - beta23*beta23
163
164 b_beta = max(0.0_rp, b_beta)
165
166 ! alpha_ij alpha_ij
167 aijaij = beta11 + beta22 + beta33
168
169 nut%x(i,1,1,e) = c*delta%x(i,1,1,e)*delta%x(i,1,1,e) &
170 * sqrt(b_beta/(aijaij + neko_eps)) &
171 * coef%mult(i,1,1,e)
172 end do
173 end do
174 !$omp end parallel do
175 if (if_corr) then
176 temperature => neko_registry%get_field_by_name(scalar_name)
177 call neko_scratch_registry%request_field(dtdx, temp_indices_buoy(1), .false.)
178 call neko_scratch_registry%request_field(dtdy, temp_indices_buoy(2), .false.)
179 call neko_scratch_registry%request_field(dtdz, temp_indices_buoy(3), .false.)
180
181 ! Calculate Richardson number
182 gmag = sqrt(vlsc2(g, g, 3))
183 if (gmag > neko_eps) then
184 n = g / gmag
185 else
186 call neko_error("The gravity vector must have at least one nonzero component")
187 endif
188 call grad(dtdx%x, dtdy%x, dtdz%x, temperature%x, coef)
189 !$omp parallel do private(e, i, j, buoyancy, du_n, du_parallel, &
190 !$omp& sh, shear_sq, ri, correction)
191 do e = 1, coef%msh%nelv
192 !OCL NORECURRENCE, NOVREC, NOALIAS
193 !DIR$ CONCURRENT
194 !GCC$ ivdep
195 do i = 1, coef%Xh%lxyz
196
197 ! Buoyancy component (numerator in Ri definition)
198 buoyancy = (-g(1) * dtdx%x(i,1,1,e) - &
199 g(2) * dtdy%x(i,1,1,e) - &
200 g(3) * dtdz%x(i,1,1,e)) / ref_temp
201
202 ! Shear component (denominator in Ri definition)
203 ! Directional derivative of velocity
204 du_n(1) = a11%x(i,1,1,e)*n(1) + a12%x(i,1,1,e)*n(2) +&
205 a13%x(i,1,1,e)*n(3)
206 du_n(2) = a21%x(i,1,1,e)*n(1) + a22%x(i,1,1,e)*n(2) +&
207 a23%x(i,1,1,e)*n(3)
208 du_n(3) = a31%x(i,1,1,e)*n(1) + a32%x(i,1,1,e)*n(2) +&
209 a33%x(i,1,1,e)*n(3)
210
211 ! Component parallel to n
212 du_parallel = du_n(1)*n(1) + du_n(2)*n(2) + du_n(3)*n(3)
213
214 ! Perpendicular (shear) components
215 do j = 1, 3
216 sh(j) = du_n(j) - du_parallel*n(j)
217 end do
218
219 ! Shear magnitude squared
220 shear_sq = sh(1)*sh(1) + sh(2)*sh(2) + sh(3)*sh(3)
221
222 ! Richardson number
223 ri = buoyancy / (shear_sq + neko_eps)
224
225 if (ri .le. ri_c) then
226 correction = sqrt(1 - ri/ri_c)
227 nut%x(i,1,1,e) = correction * nut%x(i,1,1,e)
228 else
229 nut%x(i,1,1,e) = neko_eps
230 end if
231 end do
232 end do
233 !$omp end parallel do
234 call neko_scratch_registry%relinquish_field(temp_indices_buoy)
235 end if
236
237 call coef%gs_h%op(nut, gs_op_add)
238 call col2(nut%x, coef%mult, nut%dof%size())
239
240 call neko_scratch_registry%relinquish_field(temp_indices)
241 end subroutine vreman_compute_cpu
242
243end module vreman_cpu
Compute derivative of a scalar field along a single direction.
Definition operators.f90:78
Compute the gradient of a scalar field, multiplied by the mass matrix.
Definition operators.f90:91
Coefficients.
Definition coef.f90:34
Defines a field.
Definition field.f90:34
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
Definition math.f90:60
real(kind=rp) function, public vlsc2(u, v, n)
Compute multiplication sum .
Definition math.f90:882
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:564
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:1044
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:70
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Utilities.
Definition utils.f90:35
Implements the CPU kernel for the vreman_t type.
subroutine, public vreman_compute_cpu(if_ext, t, tstep, coef, nut, delta, c, if_corr, scalar_name, ri_c, ref_temp, g)
Compute eddy viscosity on the CPU.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
field_list_t, To be able to group fields together
#define max(a, b)
Definition tensor.cu:40