Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
deardorff_cpu.f90
Go to the documentation of this file.
1! Copyright (c) 2025-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 utils, only : neko_error
36 use num_types, only : rp
37 use field_list, only : field_list_t
39 use registry, only : neko_registry
40 use field, only : field_t
41 use operators, only : strain_rate
42 use coefs, only : coef_t
43 use gs_ops, only : gs_op_add
44 use math, only : col2, neko_eps
45 use operators, only : grad
46 implicit none
47 private
48
49 public :: deardorff_compute_cpu
50
51contains
52
67 subroutine deardorff_compute_cpu(t, tstep, coef, &
68 temperature_field_name, TKE_field_name, &
69 nut, temperature_alphat, &
70 TKE_alphat, TKE_source, &
71 delta, c_k, T0, g)
72 real(kind=rp), intent(in) :: t
73 integer, intent(in) :: tstep
74 type(coef_t), intent(in) :: coef
75 character(len=*), intent(in) :: temperature_field_name
76 character(len=*), intent(in) :: tke_field_name
77 type(field_t), intent(inout) :: nut, temperature_alphat
78 type(field_t), intent(inout) :: tke_alphat, tke_source
79 type(field_t), intent(in) :: delta
80 real(kind=rp), intent(in) :: c_k, t0, g(3)
81 type(field_t), pointer :: tke, temperature
82 type(field_t), pointer :: dtdx, dtdy, dtdz
83 type(field_t), pointer :: u, v, w
84 real(kind=rp):: s11, s22, s33, s12, s13, s23
85 type(field_t), pointer :: a11, a12, a13, a21, a22, a23, a31, a32, a33
86 real(kind=rp) :: shear, buoyancy, dissipation
87 integer :: temp_indices(12)
88 real(kind=rp) :: l, n2
89 integer :: e, i
90
91 tke => neko_registry%get_field_by_name(tke_field_name)
92 temperature => neko_registry%get_field_by_name(temperature_field_name)
93
94 u => neko_registry%get_field_by_name("u")
95 v => neko_registry%get_field_by_name("v")
96 w => neko_registry%get_field_by_name("w")
97
98 call neko_scratch_registry%request_field(dtdx, temp_indices(1), .false.)
99 call neko_scratch_registry%request_field(dtdy, temp_indices(2), .false.)
100 call neko_scratch_registry%request_field(dtdz, temp_indices(3), .false.)
101
102 ! Calculate vertical temperature gradients
103 call grad(dtdx%x, dtdy%x, dtdz%x, temperature%x, coef)
104
105 call coef%gs_h%op(dtdx, gs_op_add)
106 call coef%gs_h%op(dtdy, gs_op_add)
107 call coef%gs_h%op(dtdz, gs_op_add)
108 call col2(dtdx%x, coef%mult, nut%dof%size())
109 call col2(dtdy%x, coef%mult, nut%dof%size())
110 call col2(dtdz%x, coef%mult, nut%dof%size())
111
112 ! Compute velocity gradients
113 call neko_scratch_registry%request_field(a11, temp_indices(4), .false.)
114 call neko_scratch_registry%request_field(a12, temp_indices(5), .false.)
115 call neko_scratch_registry%request_field(a13, temp_indices(6), .false.)
116 call neko_scratch_registry%request_field(a21, temp_indices(7), .false.)
117 call neko_scratch_registry%request_field(a22, temp_indices(8), .false.)
118 call neko_scratch_registry%request_field(a23, temp_indices(9), .false.)
119 call neko_scratch_registry%request_field(a31, temp_indices(10), .false.)
120 call neko_scratch_registry%request_field(a32, temp_indices(11), .false.)
121 call neko_scratch_registry%request_field(a33, temp_indices(12), .false.)
122
123 call grad(a11%x, a12%x, a13%x, u%x, coef)
124 call grad(a21%x, a22%x, a23%x, v%x, coef)
125 call grad(a31%x, a32%x, a33%x, w%x, coef)
126
127 call coef%gs_h%op(a11, gs_op_add)
128 call coef%gs_h%op(a12, gs_op_add)
129 call coef%gs_h%op(a13, gs_op_add)
130 call coef%gs_h%op(a21, gs_op_add)
131 call coef%gs_h%op(a22, gs_op_add)
132 call coef%gs_h%op(a23, gs_op_add)
133 call coef%gs_h%op(a31, gs_op_add)
134 call coef%gs_h%op(a32, gs_op_add)
135 call coef%gs_h%op(a33, gs_op_add)
136
137 call col2(a11%x, coef%mult, nut%dof%size())
138 call col2(a12%x, coef%mult, nut%dof%size())
139 call col2(a13%x, coef%mult, nut%dof%size())
140 call col2(a21%x, coef%mult, nut%dof%size())
141 call col2(a22%x, coef%mult, nut%dof%size())
142 call col2(a23%x, coef%mult, nut%dof%size())
143 call col2(a31%x, coef%mult, nut%dof%size())
144 call col2(a32%x, coef%mult, nut%dof%size())
145 call col2(a33%x, coef%mult, nut%dof%size())
146
147 ! Determine static stability and length scale
148 do concurrent(i = 1:coef%dof%size())
149 ! correct TKE if negative or nearly zero
150 if (tke%x(i,1,1,1) .lt. neko_eps) then
151 tke%x(i,1,1,1) = neko_eps
152 end if
153
154 n2 = (dtdx%x(i,1,1,1) * g(1) + &
155 dtdy%x(i,1,1,1) * g(2) + &
156 dtdz%x(i,1,1,1) * g(3)) / t0
157 if (n2 .gt. 0.0_rp) then
158 l = 0.76_rp * sqrt(tke%x(i,1,1,1) / n2)
159 l = min(l, delta%x(i,1,1,1))
160 else
161 l = delta%x(i,1,1,1)
162 end if
163
164 ! Eddy viscosity
165 nut%x(i,1,1,1) = c_k * l * sqrt(tke%x(i,1,1,1))
166
167 ! Eddy diffusivity for temperature
168 temperature_alphat%x(i,1,1,1) = (1.0_rp + 2.0_rp * l/delta%x(i,1,1,1)) &
169 * nut%x(i,1,1,1)
170 tke_alphat%x(i,1,1,1) = 2.0_rp * nut%x(i,1,1,1) ! Eddy diffusivity of TKE
171
172 s11 = a11%x(i,1,1,1) + a11%x(i,1,1,1)
173 s22 = a22%x(i,1,1,1) + a22%x(i,1,1,1)
174 s33 = a33%x(i,1,1,1) + a33%x(i,1,1,1)
175 s12 = a12%x(i,1,1,1) + a21%x(i,1,1,1)
176 s13 = a13%x(i,1,1,1) + a31%x(i,1,1,1)
177 s23 = a23%x(i,1,1,1) + a32%x(i,1,1,1)
178 ! Shear term
179 shear = nut%x(i,1,1,1) &
180 * (s11*a11%x(i,1,1,1) &
181 + s12*a12%x(i,1,1,1) &
182 + s13*a13%x(i,1,1,1) &
183 + s12*a21%x(i,1,1,1) &
184 + s22*a22%x(i,1,1,1) &
185 + s23*a23%x(i,1,1,1) &
186 + s13*a31%x(i,1,1,1) &
187 + s23*a32%x(i,1,1,1) &
188 + s33*a33%x(i,1,1,1))
189
190 ! Buoyancy term
191 buoyancy = -(g(1) * dtdx%x(i,1,1,1) + &
192 g(2) * dtdy%x(i,1,1,1) + &
193 g(3) * dtdz%x(i,1,1,1)) / t0 * temperature_alphat%x(i,1,1,1)
194
195 dissipation = -(0.19_rp + 0.74_rp * l/ delta%x(i,1,1,1)) &
196 * sqrt(tke%x(i,1,1,1)*tke%x(i,1,1,1)*tke%x(i,1,1,1)) &
197 / l
198
199 ! Add three source terms together
200 tke_source%x(i,1,1,1) = shear + buoyancy + dissipation
201 end do
202
203 call neko_scratch_registry%relinquish_field(temp_indices)
204 end subroutine deardorff_compute_cpu
205
206end module deardorff_cpu
Coefficients.
Definition coef.f90:34
Implements the CPU kernel for the deardorff_t type.
subroutine, public deardorff_compute_cpu(t, tstep, coef, temperature_field_name, tke_field_name, nut, temperature_alphat, tke_alphat, tke_source, delta, c_k, t0, g)
Compute eddy viscosity on the CPU.
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
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 strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
Compute the strain rate tensor, i.e 0.5 * du_i/dx_j + du_j/dx_i.
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
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