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
42 use json_module, only : json_file
46 use registry, only : neko_registry
47 use logger, only : log_size, neko_log
48 implicit none
49 private
50
53 type, public, extends(les_model_t) :: deardorff_t
55 real(kind=rp) :: c_k
57 real(kind=rp) :: t0
59 real(kind=rp) :: g(3)
61 character(len=:), allocatable :: temperature_field_name
63 character(len=:), allocatable :: tke_field_name
65 type(field_t), pointer :: temperature_alphat => null()
66 type(field_t), pointer :: tke_alphat => null()
68 type(field_t), pointer :: tke_source => null()
69 contains
71 procedure, pass(this) :: init => deardorff_init
73 procedure, pass(this) :: init_from_components => &
76 procedure, pass(this) :: free => deardorff_free
78 procedure, pass(this) :: compute => deardorff_compute
79 end type deardorff_t
80
81contains
85 subroutine deardorff_init(this, fluid, json)
86 class(deardorff_t), intent(inout) :: this
87 class(fluid_scheme_base_t), intent(inout), target :: fluid
88 type(json_file), intent(inout) :: json
89 character(len=:), allocatable :: temperature_field_name
90 character(len=:), allocatable :: TKE_field_name
91 character(len=:), allocatable :: nut_name
92 character(len=:), allocatable :: temperature_alphat_name, TKE_alphat_name
93 character(len=:), allocatable :: TKE_source_name
94 character(len=:), allocatable :: delta_type
95 real(kind=rp) :: c_k, t0
96 real(kind=rp), allocatable :: g(:)
97 logical :: if_ext
98 character(len=LOG_SIZE) :: log_buf
99
100 call json_get_or_default(json, "temperature_field", &
101 temperature_field_name, "temperature")
102 call json_get_or_default(json, "TKE_field", tke_field_name, "TKE")
103 call json_get_or_default(json, "nut_field", nut_name, "nut")
104 call json_get_or_default(json, "temperature_alphat_field", &
105 temperature_alphat_name, "temperature_alphat")
106 call json_get_or_default(json, "TKE_alphat_field", tke_alphat_name, &
107 "TKE_alphat")
108 call json_get_or_default(json, "TKE_source_field", tke_source_name, &
109 "TKE_source")
110 call json_get_or_default(json, "delta_type", delta_type, "pointwise")
111 call json_get_or_default(json, "c_k", c_k, 0.10_rp)
112 call json_get(json, "T0", t0)
113 call json_get_or_lookup(json, "g", g)
114 if (.not. size(g) == 3) then
115 call neko_error("The gravity vector should have 3 components")
116 end if
117 call json_get_or_default(json, "extrapolation", if_ext, .false.)
118
119 if (if_ext) then
120 call neko_error("Extrapolation is not implemented for the " // &
121 "Deardorff model.")
122 end if
123
124 call neko_log%section('LES model')
125 write(log_buf, '(A)') 'Model : deardorff'
126 call neko_log%message(log_buf)
127 write(log_buf, '(A, A)') 'Delta evaluation : ', delta_type
128 call neko_log%message(log_buf)
129 write(log_buf, '(A, E15.7)') 'c_k : ', c_k
130 call neko_log%message(log_buf)
131 write(log_buf, '(A, L1)') 'extrapolation : ', if_ext
132 call neko_log%message(log_buf)
133 call neko_log%end_section()
134
135 call deardorff_init_from_components(this, fluid, c_k, t0, &
136 temperature_field_name, tke_field_name, nut_name, &
137 temperature_alphat_name, tke_alphat_name, tke_source_name, g, &
138 delta_type, if_ext)
139
140 deallocate(temperature_field_name)
141 deallocate(tke_field_name)
142 deallocate(nut_name)
143 deallocate(temperature_alphat_name)
144 deallocate(tke_alphat_name)
145 deallocate(tke_source_name)
146 deallocate(delta_type)
147 deallocate(g)
148
149 end subroutine deardorff_init
150
165 subroutine deardorff_init_from_components(this, fluid, &
166 c_k, T0, temperature_field_name, TKE_field_name, nut_name, &
167 temperature_alphat_name, TKE_alphat_name, TKE_source_name, g, &
168 delta_type, if_ext)
169 class(deardorff_t), intent(inout) :: this
170 class(fluid_scheme_base_t), intent(inout), target :: fluid
171 real(kind=rp), intent(in) :: c_k, t0
172 character(len=*), intent(in) :: temperature_field_name
173 character(len=*), intent(in) :: TKE_field_name
174 character(len=*), intent(in) :: nut_name
175 character(len=*), intent(in) :: temperature_alphat_name, TKE_alphat_name
176 character(len=*), intent(in) :: TKE_source_name
177 real(kind=rp), intent(in) :: g(3)
178 character(len=*), intent(in) :: delta_type
179 logical, intent(in) :: if_ext
180
181 call this%free()
182 call this%init_base(fluid, nut_name, delta_type, if_ext)
183
184 call neko_registry%add_field(fluid%dm_Xh, &
185 trim(temperature_alphat_name), .true.)
186 call neko_registry%add_field(fluid%dm_Xh, trim(tke_alphat_name), .true.)
187 call neko_registry%add_field(fluid%dm_Xh, trim(tke_source_name), .true.)
188
189 this%temperature_alphat => &
190 neko_registry%get_field(trim(temperature_alphat_name))
191 this%TKE_alphat => neko_registry%get_field(trim(tke_alphat_name))
192 this%TKE_source => neko_registry%get_field(trim(tke_source_name))
193 this%c_k = c_k
194 this%T0 = t0
195 this%temperature_field_name = temperature_field_name
196 this%TKE_field_name = tke_field_name
197 this%g = -g
198
199 end subroutine deardorff_init_from_components
200
202 subroutine deardorff_free(this)
203 class(deardorff_t), intent(inout) :: this
204
205 nullify(this%temperature_alphat)
206 nullify(this%TKE_alphat)
207 nullify(this%TKE_source)
208
209 if (allocated(this%temperature_field_name)) then
210 deallocate(this%temperature_field_name)
211 end if
212
213 if (allocated(this%TKE_field_name)) then
214 deallocate(this%TKE_field_name)
215 end if
216
217 call this%free_base()
218 end subroutine deardorff_free
219
223 subroutine deardorff_compute(this, t, tstep)
224 class(deardorff_t), intent(inout) :: this
225 real(kind=rp), intent(in) :: t
226 integer, intent(in) :: tstep
227
228 type(field_t), pointer :: u, v, w, u_e, v_e, w_e
229
230 ! Compute the eddy viscosity field
231 if (neko_bcknd_device .eq. 1) then
232 call deardorff_compute_device(t, tstep, this%coef, &
233 this%temperature_field_name, this%TKE_field_name, &
234 this%nut, this%temperature_alphat, this%TKE_alphat, &
235 this%TKE_source, this%delta, this%c_k, this%T0, this%g)
236 else
237 call deardorff_compute_cpu(t, tstep, this%coef, &
238 this%temperature_field_name, this%TKE_field_name, &
239 this%nut, this%temperature_alphat, this%TKE_alphat, &
240 this%TKE_source, this%delta, this%c_k, this%T0, this%g)
241 end if
242
243 end subroutine deardorff_compute
244
245end 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:86
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:53
Base type of all fluid formulations.
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:64