Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
deardorff.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!
33!
36 use utils, only : neko_error
37 use num_types, only : rp
38 use field, only : field_t
40 use les_model, only : les_model_t
43 use json_module, only : json_file
47 use registry, only : neko_registry
48 use logger, only : log_size, neko_log
49 implicit none
50 private
51
54 type, public, extends(les_model_t) :: deardorff_t
56 real(kind=rp) :: c_k
58 real(kind=rp) :: t0
60 real(kind=rp) :: g(3)
62 character(len=:), allocatable :: temperature_field_name
64 character(len=:), allocatable :: tke_field_name
66 type(field_t), pointer :: temperature_alphat => null()
67 type(field_t), pointer :: tke_alphat => null()
69 type(field_t), pointer :: tke_source => null()
70 contains
72 procedure, pass(this) :: init => deardorff_init
74 procedure, pass(this) :: init_from_components => &
77 procedure, pass(this) :: free => deardorff_free
79 procedure, pass(this) :: compute => deardorff_compute
80 end type deardorff_t
81
82contains
86 subroutine deardorff_init(this, fluid, json)
87 class(deardorff_t), intent(inout) :: this
88 class(fluid_scheme_base_t), intent(inout), target :: fluid
89 type(json_file), intent(inout) :: json
90 character(len=:), allocatable :: temperature_field_name
91 character(len=:), allocatable :: TKE_field_name
92 character(len=:), allocatable :: nut_name
93 character(len=:), allocatable :: temperature_alphat_name, TKE_alphat_name
94 character(len=:), allocatable :: TKE_source_name
95 character(len=:), allocatable :: delta_type
96 real(kind=rp) :: c_k, t0
97 real(kind=rp), allocatable :: g(:)
98 logical :: if_ext
99 character(len=LOG_SIZE) :: log_buf
100
101 call json_get_or_default(json, "temperature_field", &
102 temperature_field_name, "temperature")
103 call json_get_or_default(json, "TKE_field", tke_field_name, "TKE")
104 call json_get_or_default(json, "nut_field", nut_name, "nut")
105 call json_get_or_default(json, "temperature_alphat_field", &
106 temperature_alphat_name, "temperature_alphat")
107 call json_get_or_default(json, "TKE_alphat_field", tke_alphat_name, &
108 "TKE_alphat")
109 call json_get_or_default(json, "TKE_source_field", tke_source_name, &
110 "TKE_source")
111 call json_get_or_default(json, "delta_type", delta_type, "pointwise")
112 call json_get_or_lookup_or_default(json, "c_k", c_k, 0.10_rp)
113 call json_get_or_lookup(json, "T0", t0)
114 call json_get_or_lookup(json, "g", g)
115 if (.not. size(g) == 3) then
116 call neko_error("The gravity vector should have 3 components")
117 end if
118 call json_get_or_default(json, "extrapolation", if_ext, .false.)
119
120 if (if_ext) then
121 call neko_error("Extrapolation is not implemented for the " // &
122 "Deardorff model.")
123 end if
124
125 call neko_log%section('LES model')
126 write(log_buf, '(A)') 'Model : deardorff'
127 call neko_log%message(log_buf)
128 write(log_buf, '(A, A)') 'Delta evaluation : ', delta_type
129 call neko_log%message(log_buf)
130 write(log_buf, '(A, E15.7)') 'c_k : ', c_k
131 call neko_log%message(log_buf)
132 write(log_buf, '(A, L1)') 'extrapolation : ', if_ext
133 call neko_log%message(log_buf)
134 call neko_log%end_section()
135
136 call deardorff_init_from_components(this, fluid, c_k, t0, &
137 temperature_field_name, tke_field_name, nut_name, &
138 temperature_alphat_name, tke_alphat_name, tke_source_name, g, &
139 delta_type, if_ext)
140
141 deallocate(temperature_field_name)
142 deallocate(tke_field_name)
143 deallocate(nut_name)
144 deallocate(temperature_alphat_name)
145 deallocate(tke_alphat_name)
146 deallocate(tke_source_name)
147 deallocate(delta_type)
148 deallocate(g)
149
150 end subroutine deardorff_init
151
166 subroutine deardorff_init_from_components(this, fluid, &
167 c_k, T0, temperature_field_name, TKE_field_name, nut_name, &
168 temperature_alphat_name, TKE_alphat_name, TKE_source_name, g, &
169 delta_type, if_ext)
170 class(deardorff_t), intent(inout) :: this
171 class(fluid_scheme_base_t), intent(inout), target :: fluid
172 real(kind=rp), intent(in) :: c_k, t0
173 character(len=*), intent(in) :: temperature_field_name
174 character(len=*), intent(in) :: TKE_field_name
175 character(len=*), intent(in) :: nut_name
176 character(len=*), intent(in) :: temperature_alphat_name, TKE_alphat_name
177 character(len=*), intent(in) :: TKE_source_name
178 real(kind=rp), intent(in) :: g(3)
179 character(len=*), intent(in) :: delta_type
180 logical, intent(in) :: if_ext
181
182 call this%free()
183 call this%init_base(fluid, nut_name, delta_type, if_ext)
184
185 call neko_registry%add_field(fluid%dm_Xh, &
186 trim(temperature_alphat_name), .true.)
187 call neko_registry%add_field(fluid%dm_Xh, trim(tke_alphat_name), .true.)
188 call neko_registry%add_field(fluid%dm_Xh, trim(tke_source_name), .true.)
189
190 this%temperature_alphat => &
191 neko_registry%get_field(trim(temperature_alphat_name))
192 this%TKE_alphat => neko_registry%get_field(trim(tke_alphat_name))
193 this%TKE_source => neko_registry%get_field(trim(tke_source_name))
194 this%c_k = c_k
195 this%T0 = t0
196 this%temperature_field_name = temperature_field_name
197 this%TKE_field_name = tke_field_name
198 this%g = -g
199
200 end subroutine deardorff_init_from_components
201
203 subroutine deardorff_free(this)
204 class(deardorff_t), intent(inout) :: this
205
206 nullify(this%temperature_alphat)
207 nullify(this%TKE_alphat)
208 nullify(this%TKE_source)
209
210 if (allocated(this%temperature_field_name)) then
211 deallocate(this%temperature_field_name)
212 end if
213
214 if (allocated(this%TKE_field_name)) then
215 deallocate(this%TKE_field_name)
216 end if
217
218 call this%free_base()
219 end subroutine deardorff_free
220
224 subroutine deardorff_compute(this, t, tstep)
225 class(deardorff_t), intent(inout) :: this
226 real(kind=rp), intent(in) :: t
227 integer, intent(in) :: tstep
228
229 type(field_t), pointer :: u, v, w, u_e, v_e, w_e
230
231 ! Compute the eddy viscosity field
232 if (neko_bcknd_device .eq. 1) then
233 call deardorff_compute_device(t, tstep, this%coef, &
234 this%temperature_field_name, this%TKE_field_name, &
235 this%nut, this%temperature_alphat, this%TKE_alphat, &
236 this%TKE_source, this%delta, this%c_k, this%T0, this%g)
237 else
238 call deardorff_compute_cpu(t, tstep, this%coef, &
239 this%temperature_field_name, this%TKE_field_name, &
240 this%nut, this%temperature_alphat, this%TKE_alphat, &
241 this%TKE_source, this%delta, this%c_k, this%T0, this%g)
242 end if
243
244 end subroutine deardorff_compute
245
246end module deardorff
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
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.
Implements the device kernel for the deardorff_t type.
subroutine, public deardorff_compute_device(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 device.
Implements deardorff_t.
Definition deardorff.f90:35
subroutine deardorff_init(this, fluid, json)
Constructor.
Definition deardorff.f90:87
subroutine deardorff_free(this)
Destructor for the deardorff_t class.
subroutine deardorff_compute(this, t, tstep)
Compute eddy viscosity.
subroutine deardorff_init_from_components(this, fluid, c_k, t0, temperature_field_name, tke_field_name, nut_name, temperature_alphat_name, tke_alphat_name, tke_source_name, g, delta_type, if_ext)
Constructor from components.
Defines a field.
Definition field.f90:34
Utilities for retrieving parameters from the case files.
Implements les_model_t.
Definition les_model.f90:35
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
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
Utilities.
Definition utils.f90:35
Implements the deardorff LES model.
Definition deardorff.f90:54
Base type of all fluid formulations.
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:64