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 do concurrent(e = 1:coef%msh%nelv)
141 do concurrent(i = 1:coef%Xh%lxyz)
142 ! beta_ij = alpha_mi alpha_mj
143 beta11 = a11%x(i,1,1,e)**2 + a21%x(i,1,1,e)**2 + a31%x(i,1,1,e)**2
144 beta22 = a12%x(i,1,1,e)**2 + a22%x(i,1,1,e)**2 + a32%x(i,1,1,e)**2
145 beta33 = a13%x(i,1,1,e)**2 + a23%x(i,1,1,e)**2 + a33%x(i,1,1,e)**2
146 beta12 = a11%x(i,1,1,e)*a12%x(i,1,1,e) + &
147 a21%x(i,1,1,e)*a22%x(i,1,1,e) + &
148 a31%x(i,1,1,e)*a32%x(i,1,1,e)
149 beta13 = a11%x(i,1,1,e)*a13%x(i,1,1,e) + &
150 a21%x(i,1,1,e)*a23%x(i,1,1,e) + &
151 a31%x(i,1,1,e)*a33%x(i,1,1,e)
152 beta23 = a12%x(i,1,1,e)*a13%x(i,1,1,e) + &
153 a22%x(i,1,1,e)*a23%x(i,1,1,e) + &
154 a32%x(i,1,1,e)*a33%x(i,1,1,e)
155
156 b_beta = beta11*beta22 - beta12*beta12 + beta11*beta33 &
157 - beta13*beta13 + beta22*beta33 - beta23*beta23
158
159 b_beta = max(0.0_rp, b_beta)
160
161 ! alpha_ij alpha_ij
162 aijaij = beta11 + beta22 + beta33
163
164 nut%x(i,1,1,e) = c*delta%x(i,1,1,e)*delta%x(i,1,1,e) &
165 * sqrt(b_beta/(aijaij + neko_eps)) &
166 * coef%mult(i,1,1,e)
167 end do
168 end do
169 if (if_corr) then
170 temperature => neko_registry%get_field_by_name(scalar_name)
171 call neko_scratch_registry%request_field(dtdx, temp_indices_buoy(1), .false.)
172 call neko_scratch_registry%request_field(dtdy, temp_indices_buoy(2), .false.)
173 call neko_scratch_registry%request_field(dtdz, temp_indices_buoy(3), .false.)
174
175 ! Calculate Richardson number
176 gmag = sqrt(vlsc2(g, g, 3))
177 if (gmag > neko_eps) then
178 n = g / gmag
179 else
180 call neko_error("The gravity vector must have at least one nonzero component")
181 endif
182 call grad(dtdx%x, dtdy%x, dtdz%x, temperature%x, coef)
183 do concurrent(e = 1:coef%msh%nelv)
184 do concurrent(i = 1:coef%Xh%lxyz)
185
186 ! Buoyancy component (numerator in Ri definition)
187 buoyancy = (-g(1) * dtdx%x(i,1,1,e) - &
188 g(2) * dtdy%x(i,1,1,e) - &
189 g(3) * dtdz%x(i,1,1,e)) / ref_temp
190
191 ! Shear component (denominator in Ri definition)
192 ! Directional derivative of velocity
193 du_n(1) = a11%x(i,1,1,e)*n(1) + a12%x(i,1,1,e)*n(2) +&
194 a13%x(i,1,1,e)*n(3)
195 du_n(2) = a21%x(i,1,1,e)*n(1) + a22%x(i,1,1,e)*n(2) +&
196 a23%x(i,1,1,e)*n(3)
197 du_n(3) = a31%x(i,1,1,e)*n(1) + a32%x(i,1,1,e)*n(2) +&
198 a33%x(i,1,1,e)*n(3)
199
200 ! Component parallel to n
201 du_parallel = du_n(1)*n(1) + du_n(2)*n(2) + du_n(3)*n(3)
202
203 ! Perpendicular (shear) components
204 do concurrent(j = 1:3)
205 sh(j) = du_n(j) - du_parallel*n(j)
206 end do
207
208 ! Shear magnitude squared
209 shear_sq = sh(1)*sh(1) + sh(2)*sh(2) + sh(3)*sh(3)
210
211 ! Richardson number
212 ri = buoyancy / (shear_sq + neko_eps)
213
214 if (ri .le. ri_c) then
215 correction = sqrt(1 - ri/ri_c)
216 nut%x(i,1,1,e) = correction * nut%x(i,1,1,e)
217 else
218 nut%x(i,1,1,e) = neko_eps
219 end if
220 end do
221 end do
222 call neko_scratch_registry%relinquish_field(temp_indices_buoy)
223 end if
224
225 call coef%gs_h%op(nut, gs_op_add)
226 call col2(nut%x, coef%mult, nut%dof%size())
227
228 call neko_scratch_registry%relinquish_field(temp_indices)
229 end subroutine vreman_compute_cpu
230
231end module vreman_cpu
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:712
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:462
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:69
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public grad(ux, uy, uz, u, coef)
Compute the gradient of a scalar field.
subroutine, public dudxyz(du, u, dr, ds, dt, coef)
Compute derivative of a scalar field along a single direction.
Definition operators.f90:92
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:58
field_list_t, To be able to group fields together
#define max(a, b)
Definition tensor.cu:40