Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
les_model.f90
Go to the documentation of this file.
1! Copyright (c) 2023-2024, 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 num_types, only : rp
38 use fluid_pnpn, only : fluid_pnpn_t
40 use rhs_maker, only : rhs_maker_sumab_t, rhs_maker_sumab_fctry
41 use field, only : field_t
43 use json_module, only : json_file
45 use dofmap, only : dofmap_t
46 use coefs, only : coef_t
47 use gs_ops, only : gs_op_add
50 use math, only : col2
51 use device_math, only : device_col2
54 use comm, only : pe_rank
55 implicit none
56 private
57
58 ! List of all possible types created by the factory routine
59 character(len=20) :: delta_known_types(3) = [character(len=20) :: &
60 "pointwise", &
61 "elementwise_average", &
62 "elementwise_max"]
63
65 type, abstract, public :: les_model_t
67 type(time_scheme_controller_t), pointer :: ext_bdf => null()
69 type(field_series_t), pointer :: ulag => null()
70 type(field_series_t), pointer :: vlag => null()
71 type(field_series_t), pointer :: wlag => null()
73 class(rhs_maker_sumab_t), allocatable :: sumab
75 logical :: if_ext = .false.
77 type(field_t), pointer :: nut => null()
79 character(len=:), allocatable :: delta_type
81 type(field_t), pointer :: delta => null()
83 type(coef_t), pointer :: coef => null()
84 contains
86 procedure, pass(this) :: init_base => les_model_init_base
88 procedure, pass(this) :: free_base => les_model_free_base
90 procedure, pass(this) :: compute_delta => les_model_compute_delta
92 procedure(les_model_init), pass(this), deferred :: init
94 procedure(les_model_free), pass(this), deferred :: free
96 procedure(les_model_compute), pass(this), deferred :: compute
97 end type les_model_t
98
99 abstract interface
100
103 subroutine les_model_compute(this, t, tstep)
104 import les_model_t, rp
105 class(les_model_t), intent(inout) :: this
106 real(kind=rp), intent(in) :: t
107 integer, intent(in) :: tstep
108 end subroutine les_model_compute
109 end interface
110
111 abstract interface
112
115 subroutine les_model_init(this, fluid, json)
116 import les_model_t, json_file, fluid_scheme_base_t
117 class(les_model_t), intent(inout) :: this
118 class(fluid_scheme_base_t), intent(inout), target :: fluid
119 type(json_file), intent(inout) :: json
120 end subroutine les_model_init
121 end interface
122
123 abstract interface
124
125 subroutine les_model_free(this)
126 import les_model_t
127 class(les_model_t), intent(inout) :: this
128 end subroutine les_model_free
129 end interface
130
131 interface
132
135 module subroutine les_model_allocator(object, type_name)
136 class(les_model_t), allocatable, intent(inout) :: object
137 character(len=*), intent(in) :: type_name
138 end subroutine les_model_allocator
139 end interface
140
141 interface
142
149 module subroutine les_model_factory(object, type_name, fluid, json)
150 class(les_model_t), allocatable, intent(inout) :: object
151 character(len=*), intent(in) :: type_name
152 class(fluid_scheme_base_t), intent(inout) :: fluid
153 type(json_file), intent(inout) :: json
154 end subroutine les_model_factory
155 end interface
156
157 !
158 ! Machinery for injecting user-defined types
159 !
160
164 abstract interface
165 subroutine les_model_allocate(obj)
166 import les_model_t
167 class(les_model_t), allocatable, intent(inout) :: obj
168 end subroutine les_model_allocate
169 end interface
170
171 interface
172
173 module subroutine register_les_model(type_name, allocator)
174 character(len=*), intent(in) :: type_name
175 procedure(les_model_allocate), pointer, intent(in) :: allocator
176 end subroutine register_les_model
177 end interface
178
179 ! A name-allocator pair for user-defined types. A helper type to define a
180 ! registry of custom allocators.
181 type allocator_entry
182 character(len=20) :: type_name
183 procedure(les_model_allocate), pointer, nopass :: allocator
184 end type allocator_entry
185
187 type(allocator_entry), allocatable :: les_model_registry(:)
188
190 integer :: les_model_registry_size = 0
191
192 public :: les_model_factory, les_model_allocator, register_les_model, &
193 les_model_allocate
194
195
196contains
202 subroutine les_model_init_base(this, fluid, nut_name, delta_type, if_ext)
203 class(les_model_t), intent(inout) :: this
204 class(fluid_scheme_base_t), intent(inout), target :: fluid
205 character(len=*), intent(in) :: nut_name
206 character(len=*), intent(in) :: delta_type
207 logical, intent(in) :: if_ext
208
209 associate(dofmap => fluid%dm_Xh, &
210 coef => fluid%c_Xh)
211
212 call neko_field_registry%add_field(dofmap, trim(nut_name), .true.)
213 call neko_field_registry%add_field(dofmap, "les_delta", .true.)
214 this%nut => neko_field_registry%get_field(trim(nut_name))
215 this%delta => neko_field_registry%get_field("les_delta")
216 this%coef => coef
217 this%delta_type = delta_type
218 this%if_ext = if_ext
219
220 if (pe_rank .eq. 0) then
221 if (if_ext .eqv. .true.) then
222 call neko_warning("Extrapolation of the velocity in eddy &
223 &viscosity estimation might be unstable.")
224 else
225 call neko_warning("The time integration for eddy viscosity &
226 &estimation is only first-order accurate")
227 end if
228 end if
229
230 call this%compute_delta()
231
232 select type (fluid)
233 type is (fluid_pnpn_t)
234 this%ulag => fluid%ulag
235 this%vlag => fluid%vlag
236 this%wlag => fluid%wlag
237 ! Setup backend dependent summation of AB/BDF
238 this%ext_bdf => fluid%ext_bdf
239 call rhs_maker_sumab_fctry(this%sumab)
240 class default
241 if (this%if_ext .eqv. .true.) then
242 call neko_error("Fluid scheme does not support &
243 &velocity extrapolation")
244 end if
245 end select
246
247 end associate
248 end subroutine les_model_init_base
249
251 subroutine les_model_free_base(this)
252 class(les_model_t), intent(inout) :: this
253
254 nullify(this%nut)
255 nullify(this%delta)
256 nullify(this%coef)
257 if (allocated(this%sumab)) then
258 deallocate(this%sumab)
259 end if
260 end subroutine les_model_free_base
261
268 subroutine les_model_compute_delta(this)
269 class(les_model_t), intent(inout) :: this
270 integer :: e, i, j, k
271 integer :: im, ip, jm, jp, km, kp
272 real(kind=rp) :: di, dj, dk, ndim_inv, volume_element
273 integer :: lx_half, ly_half, lz_half
274 character(len=:), allocatable :: type_string
275
276 lx_half = this%coef%Xh%lx / 2
277 ly_half = this%coef%Xh%ly / 2
278 lz_half = this%coef%Xh%lz / 2
279
280 if (this%delta_type .eq. "elementwise_max") then
281 ! use a same length scale throughout an entire element
282 ! the length scale is based on maximum GLL spacing
283 do e = 1, this%coef%msh%nelv
284 di = (this%coef%dof%x(lx_half, 1, 1, e) &
285 - this%coef%dof%x(lx_half + 1, 1, 1, e))**2 &
286 + (this%coef%dof%y(lx_half, 1, 1, e) &
287 - this%coef%dof%y(lx_half + 1, 1, 1, e))**2 &
288 + (this%coef%dof%z(lx_half, 1, 1, e) &
289 - this%coef%dof%z(lx_half + 1, 1, 1, e))**2
290
291 dj = (this%coef%dof%x(1, ly_half, 1, e) &
292 - this%coef%dof%x(1, ly_half + 1, 1, e))**2 &
293 + (this%coef%dof%y(1, ly_half, 1, e) &
294 - this%coef%dof%y(1, ly_half + 1, 1, e))**2 &
295 + (this%coef%dof%z(1, ly_half, 1, e) &
296 - this%coef%dof%z(1, ly_half + 1, 1, e))**2
297
298 dk = (this%coef%dof%x(1, 1, lz_half, e) &
299 - this%coef%dof%x(1, 1, lz_half + 1, e))**2 &
300 + (this%coef%dof%y(1, 1, lz_half, e) &
301 - this%coef%dof%y(1, 1, lz_half + 1, e))**2 &
302 + (this%coef%dof%z(1, 1, lz_half, e) &
303 - this%coef%dof%z(1, 1, lz_half + 1, e))**2
304 di = sqrt(di)
305 dj = sqrt(dj)
306 dk = sqrt(dk)
307 this%delta%x(:,:,:,e) = (di * dj * dk)**(1.0_rp / 3.0_rp)
308 end do
309 else if (this%delta_type .eq. "elementwise_average") then
310 ! use a same length scale throughout an entire element
311 ! the length scale is based on (volume)^(1/3)/(N+1)
312 do e = 1, this%coef%msh%nelv
313 volume_element = 0.0_rp
314 do k = 1, this%coef%Xh%lx * this%coef%Xh%ly * this%coef%Xh%lz
315 volume_element = volume_element + this%coef%B(k, 1, 1, e)
316 end do
317 this%delta%x(:,:,:,e) = (volume_element / this%coef%Xh%lx &
318 / this%coef%Xh%ly / this%coef%Xh%lz)**(1.0_rp / 3.0_rp)
319 end do
320 else if (this%delta_type .eq. "pointwise") then
321 do e = 1, this%coef%msh%nelv
322 do k = 1, this%coef%Xh%lz
323 km = max(1, k-1)
324 kp = min(this%coef%Xh%lz, k+1)
325
326 do j = 1, this%coef%Xh%ly
327 jm = max(1, j-1)
328 jp = min(this%coef%Xh%ly, j+1)
329
330 do i = 1, this%coef%Xh%lx
331 im = max(1, i-1)
332 ip = min(this%coef%Xh%lx, i+1)
333
334 di = (this%coef%dof%x(ip, j, k, e) - &
335 this%coef%dof%x(im, j, k, e))**2 &
336 + (this%coef%dof%y(ip, j, k, e) - &
337 this%coef%dof%y(im, j, k, e))**2 &
338 + (this%coef%dof%z(ip, j, k, e) - &
339 this%coef%dof%z(im, j, k, e))**2
340
341 dj = (this%coef%dof%x(i, jp, k, e) - &
342 this%coef%dof%x(i, jm, k, e))**2 &
343 + (this%coef%dof%y(i, jp, k, e) - &
344 this%coef%dof%y(i, jm, k, e))**2 &
345 + (this%coef%dof%z(i, jp, k, e) - &
346 this%coef%dof%z(i, jm, k, e))**2
347
348 dk = (this%coef%dof%x(i, j, kp, e) - &
349 this%coef%dof%x(i, j, km, e))**2 &
350 + (this%coef%dof%y(i, j, kp, e) - &
351 this%coef%dof%y(i, j, km, e))**2 &
352 + (this%coef%dof%z(i, j, kp, e) - &
353 this%coef%dof%z(i, j, km, e))**2
354
355 di = sqrt(di) / (ip - im)
356 dj = sqrt(dj) / (jp - jm)
357 dk = sqrt(dk) / (kp - km)
358 this%delta%x(i,j,k,e) = (di * dj * dk)**(1.0_rp / 3.0_rp)
359
360 end do
361 end do
362 end do
363 end do
364 else
365 call neko_type_error("delta_type for LES model", &
366 this%delta_type, delta_known_types)
367 stop
368 end if
369
370 if (neko_bcknd_device .eq. 1) then
371 call device_memcpy(this%delta%x, this%delta%x_d, this%delta%dof%size(),&
372 host_to_device, sync = .false.)
373 call this%coef%gs_h%op(this%delta%x, this%delta%dof%size(), gs_op_add)
374 call device_col2(this%delta%x_d, this%coef%mult_d, this%delta%dof%size())
375 else
376 call this%coef%gs_h%op(this%delta%x, this%delta%dof%size(), gs_op_add)
377 call col2(this%delta%x, this%coef%mult, this%delta%dof%size())
378 end if
379
380 end subroutine les_model_compute_delta
381
382end module les_model
Copy data between host and device (or device and device)
Definition device.F90:65
Compute eddy viscosity.
Common constructor.
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
integer pe_rank
MPI rank.
Definition comm.F90:51
subroutine, public device_col2(a_d, b_d, n)
Vector multiplication .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:46
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Stores a series fields.
Defines a field.
Definition field.f90:34
Modular version of the Classic Nek5000 Pn/Pn formulation for fluids.
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
Implements les_model_t.
Definition les_model.f90:35
character(len=20), dimension(3) delta_known_types
Definition les_model.f90:59
subroutine les_model_compute_delta(this)
Compute the LES lengthscale. For each GLL point, we take the distance between its neighbours in all 3...
subroutine les_model_init_base(this, fluid, nut_name, delta_type, if_ext)
Constructor for the les_model_t (base) class.
subroutine les_model_free_base(this)
Destructor for the les_model_t (base) class.
Definition math.f90:60
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:728
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
Definition rhs_maker.f90:38
Compound scheme for the advection and diffusion operators in a transport equation.
Utilities.
Definition utils.f90:35
character(:) function, allocatable, public concat_string_array(array, sep, prepend)
Concatenate an array of strings into one string with array items separated by spaces.
Definition utils.f90:276
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:266
subroutine, public neko_type_error(base_type, wrong_type, known_types)
Reports an error allocating a type for a particular base pointer class.
Definition utils.f90:250
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Base type of all fluid formulations.
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:65
Abstract type to compute extrapolated velocity field for the pressure equation.
Definition rhs_maker.f90:46
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
#define max(a, b)
Definition tensor.cu:40