Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
vreman.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!
33!
35module vreman
36 use num_types, only : rp
37 use les_model, only : les_model_t
38 use field, only : field_t
42 use json_module, only : json_file
44 use utils, only : neko_error
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) :: vreman_t
56 real(kind=rp) :: c
57 real(kind=rp) :: ri_c, ref_temp
58 real(kind=rp) :: g(3)
59 logical :: if_corr
60 character(len=:), allocatable :: scalar_name
61 contains
63 procedure, pass(this) :: init => vreman_init
65 procedure, pass(this) :: init_from_components => &
68 procedure, pass(this) :: free => vreman_free
70 procedure, pass(this) :: compute => vreman_compute
71 end type vreman_t
72
73contains
77 subroutine vreman_init(this, fluid, json)
78 class(vreman_t), intent(inout) :: this
79 class(fluid_scheme_base_t), intent(inout), target :: fluid
80 type(json_file), intent(inout) :: json
81 character(len=:), allocatable :: nut_name
82
83 real(kind=rp) :: c, ri_c, ref_temp
84 real(kind=rp), allocatable :: g(:)
85 character(len=:), allocatable :: scalar_name
86 character(len=:), allocatable :: delta_type
87 logical :: if_ext, if_corr
88 character(len=LOG_SIZE) :: log_buf
89
90 call json_get_or_default(json, "nut_field", nut_name, "nut")
91 call json_get_or_default(json, "delta_type", delta_type, "pointwise")
92 ! Based on the Smagorinsky Cs = 0.17.
93 call json_get_or_lookup_or_default(json, "c", c, 0.07_rp)
94 call json_get_or_default(json, "extrapolation", if_ext, .false.)
95 call json_get_or_default(json, "buoyancy_correction", if_corr, .false.)
96
97 if (if_corr) then
98 call json_get(json, "scalar_field", scalar_name)
99 call json_get_or_lookup(json, "Ri_c", ri_c)
100 call json_get_or_lookup(json, "reference_temperature", ref_temp)
101 call json_get_or_lookup(json, "g", g)
102 if (.not. size(g) == 3) then
103 call neko_error("The gravity vector should have 3 components")
104 end if
105
106 if (ri_c .le. 0.0_rp) then
107 call neko_error("The critical Richardson number should be positive.")
108 end if
109
110 if (ref_temp .le. 0.0_rp) then
111 call neko_error("The reference temperature should be positive.")
112 end if
113 else
114 ! These values are not used
115 g = [0.0_rp, 0.0_rp, 0.0_rp]
116 ri_c = 0.0_rp
117 ref_temp = 0.0_rp
118 scalar_name = ""
119 end if
120
121 call neko_log%section('LES model')
122 write(log_buf, '(A)') 'Model : Vreman'
123 call neko_log%message(log_buf)
124 write(log_buf, '(A, A)') 'Delta evaluation : ', delta_type
125 call neko_log%message(log_buf)
126 write(log_buf, '(A, E15.7)') 'c : ', c
127 call neko_log%message(log_buf)
128 write(log_buf, '(A, L1)') 'extrapolation : ', if_ext
129 call neko_log%message(log_buf)
130 write(log_buf, '(A, L1)') 'buoyancy correction : ', if_corr
131 call neko_log%message(log_buf)
132 call neko_log%end_section()
133
134 call vreman_init_from_components(this, fluid, c, nut_name, &
135 delta_type, if_ext, if_corr, scalar_name, ri_c, ref_temp, g)
136
137 end subroutine vreman_init
138
150 subroutine vreman_init_from_components(this, fluid, c, nut_name, &
151 delta_type, if_ext, if_corr, scalar_name, ri_c, ref_temp, g)
152 class(vreman_t), intent(inout) :: this
153 class(fluid_scheme_base_t), intent(inout), target :: fluid
154 real(kind=rp) :: c, ri_c, ref_temp
155 real(kind=rp) :: g(3)
156 character(len=*), intent(in) :: scalar_name
157 character(len=*), intent(in) :: nut_name
158 character(len=*), intent(in) :: delta_type
159 logical, intent(in) :: if_ext, if_corr
160
161 call this%free()
162
163 call this%init_base(fluid, nut_name, delta_type, if_ext)
164 this%c = c
165 this%if_corr = if_corr
166 this%ri_c = ri_c
167 this%ref_temp = ref_temp
168 this%g = g
169 this%scalar_name = scalar_name
170
171 end subroutine vreman_init_from_components
172
174 subroutine vreman_free(this)
175 class(vreman_t), intent(inout) :: this
176
177 call this%free_base()
178 end subroutine vreman_free
179
183 subroutine vreman_compute(this, t, tstep)
184 class(vreman_t), intent(inout) :: this
185 real(kind=rp), intent(in) :: t
186 integer, intent(in) :: tstep
187
188 type(field_t), pointer :: u, v, w, u_e, v_e, w_e
189
190 if (this%if_ext) then
191 ! Extrapolate the velocity fields
192 associate(ulag => this%ulag, vlag => this%vlag, &
193 wlag => this%wlag, ext_bdf => this%ext_bdf)
194
195 u => neko_registry%get_field_by_name("u")
196 v => neko_registry%get_field_by_name("v")
197 w => neko_registry%get_field_by_name("w")
198 u_e => neko_registry%get_field_by_name("u_e")
199 v_e => neko_registry%get_field_by_name("v_e")
200 w_e => neko_registry%get_field_by_name("w_e")
201
202 call this%sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
203 ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
204
205 end associate
206 end if
207
208 ! Compute the eddy viscosity field
209 if (neko_bcknd_device .eq. 1) then
210 call vreman_compute_device(this%if_ext, t, tstep, this%coef, &
211 this%nut, this%delta, this%c, this%if_corr, &
212 this%scalar_name, this%ri_c, this%ref_temp, this%g)
213 else
214 call vreman_compute_cpu(this%if_ext, t, tstep, this%coef, &
215 this%nut, this%delta, this%c, this%if_corr, &
216 this%scalar_name, this%ri_c, this%ref_temp, this%g)
217 end if
218
219 end subroutine vreman_compute
220
221end module vreman
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.
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 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.
Implements the device kernel for the vreman_t type.
subroutine, public vreman_compute_device(if_ext, t, tstep, coef, nut, delta, c, if_corr, scalar_name, ri_c, ref_temp, g)
Compute eddy viscosity on the device.
Implements vreman_t.
Definition vreman.f90:35
subroutine vreman_free(this)
Destructor for the les_model_t (base) class.
Definition vreman.f90:175
subroutine vreman_init_from_components(this, fluid, c, nut_name, delta_type, if_ext, if_corr, scalar_name, ri_c, ref_temp, g)
Constructor from components.
Definition vreman.f90:152
subroutine vreman_compute(this, t, tstep)
Compute eddy viscosity.
Definition vreman.f90:184
subroutine vreman_init(this, fluid, json)
Constructor.
Definition vreman.f90:78
Base type of all fluid formulations.
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:64
Implements the Vreman LES model.
Definition vreman.f90:54