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
dynamic_smagorinsky.f90
Go to the documentation of this file.
1! Copyright (c) 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!
35 use num_types, only : rp
36 use field, only : field_t
38 use les_model, only : les_model_t
40 use json_module, only : json_file
41 use utils, only : neko_error
45 use logger, only : log_size, neko_log
48 implicit none
49 private
50
53 type, public, extends(les_model_t) :: dynamic_smagorinsky_t
55 type(field_t) :: c_dyn
57 type(elementwise_filter_t) :: test_filter
61 type(field_t) :: mij(6)
63 type(field_t) :: lij(6)
65 type(field_t) :: num
67 type(field_t) :: den
68 contains
70 procedure, pass(this) :: init => dynamic_smagorinsky_init
72 procedure, pass(this) :: free => dynamic_smagorinsky_free
74 procedure, pass(this) :: compute => dynamic_smagorinsky_compute
75
77
78contains
82 subroutine dynamic_smagorinsky_init(this, fluid, json)
83 class(dynamic_smagorinsky_t), intent(inout) :: this
84 class(fluid_scheme_base_t), intent(inout), target :: fluid
85 type(json_file), intent(inout) :: json
86 character(len=:), allocatable :: nut_name
87 integer :: i
88 character(len=:), allocatable :: delta_type
89 logical :: if_ext
90 character(len=LOG_SIZE) :: log_buf
91
92 associate(dofmap => fluid%dm_Xh, &
93 coef => fluid%c_Xh)
94
95 call json_get_or_default(json, "nut_field", nut_name, "nut")
96 call json_get_or_default(json, "delta_type", delta_type, "pointwise")
97 call json_get_or_default(json, "extrapolation", if_ext, .true.)
98
99 call this%free()
100 call this%init_base(fluid, nut_name, delta_type, if_ext)
101 call this%test_filter%init(json, coef)
102 call set_ds_filt(this%test_filter)
103
104 call neko_log%section('LES model')
105 write(log_buf, '(A)') 'Model : Dynamic Smagorinsky'
106 call neko_log%message(log_buf)
107 write(log_buf, '(A, A)') 'Delta evaluation : ', delta_type
108 call neko_log%message(log_buf)
109 write(log_buf, '(A, A)') 'Test filter type : ', &
110 this%test_filter%filter_type
111 call neko_log%message(log_buf)
112 write(log_buf, '(A, L1)') 'extrapolation : ', if_ext
113 call neko_log%message(log_buf)
114 call neko_log%end_section()
115
116 call this%c_dyn%init(dofmap, "ds_c_dyn")
117 call this%num%init(dofmap, "ds_num")
118 call this%den%init(dofmap, "ds_den")
119
120 do i = 1, 6
121 call this%mij(i)%init(dofmap)
122 call this%lij(i)%init(dofmap)
123 end do
124
125 end associate
126
127 end subroutine dynamic_smagorinsky_init
128
131 class(dynamic_smagorinsky_t), intent(inout) :: this
132 integer :: i
133
134 call this%c_dyn%free()
135 do i = 1, 6
136 call this%mij(i)%free()
137 call this%lij(i)%free()
138 end do
139 call this%num%free()
140 call this%den%free()
141 call this%test_filter%free()
142 call this%free_base()
143
144 end subroutine dynamic_smagorinsky_free
145
149 subroutine dynamic_smagorinsky_compute(this, t, tstep)
150 class(dynamic_smagorinsky_t), intent(inout) :: this
151 real(kind=rp), intent(in) :: t
152 integer, intent(in) :: tstep
153
154 type(field_t), pointer :: u, v, w, u_e, v_e, w_e
155
156 if (this%if_ext .eqv. .true.) then
157 ! Extrapolate the velocity fields
158 associate(ulag => this%ulag, vlag => this%vlag, &
159 wlag => this%wlag, ext_bdf => this%ext_bdf)
160
161 u => neko_field_registry%get_field_by_name("u")
162 v => neko_field_registry%get_field_by_name("v")
163 w => neko_field_registry%get_field_by_name("w")
164 u_e => neko_field_registry%get_field_by_name("u_e")
165 v_e => neko_field_registry%get_field_by_name("v_e")
166 w_e => neko_field_registry%get_field_by_name("w_e")
167
168 call this%sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
169 ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
170
171 end associate
172 end if
173
174 ! Compute the eddy viscosity field
175 if (neko_bcknd_device .eq. 1) then
176 call dynamic_smagorinsky_compute_device(this%if_ext, t, tstep, &
177 this%coef, this%nut, &
178 this%delta, this%c_dyn, this%test_filter, &
179 this%mij, this%lij, this%num, this%den)
180 else
181 call dynamic_smagorinsky_compute_cpu(this%if_ext, t, tstep, &
182 this%coef, this%nut, &
183 this%delta, this%c_dyn, this%test_filter, &
184 this%mij, this%lij, this%num, this%den)
185 end if
186
187 end subroutine dynamic_smagorinsky_compute
188
190 subroutine set_ds_filt(filter_1d)
191 type(elementwise_filter_t), intent(inout) :: filter_1d
192 integer :: i
193
194 if (filter_1d%nx .le. 2) then
195 call neko_error("Dynamic Smagorinsky model error: test filter is not &
196 &defined for the current polynomial order")
197 end if
198 if (mod(filter_1d%nx,2) .eq. 0) then ! number of grid spacing is odd
199 ! cutoff at polynomial order int((filter_1d%nx)/2)
200 filter_1d%trnsfr(int((filter_1d%nx)/2)) = 0.95_rp
201 filter_1d%trnsfr(int((filter_1d%nx)/2) + 1) = 0.50_rp
202 filter_1d%trnsfr(int((filter_1d%nx)/2) + 2) = 0.05_rp
203 if ((int((filter_1d%nx)/2) + 2) .lt. filter_1d%nx) then
204 do i = int((filter_1d%nx)/2) + 3, filter_1d%nx
205 filter_1d%trnsfr(i) = 0.0_rp
206 end do
207 end if
208 ! make delta_ratio = (nx-1)/(nt-1) as close to 2
209 filter_1d%nt = int(filter_1d%nx/2) + 1
210 else ! number of grid spacing is even
211 ! cutoff at polynomial order int((filter_1d%nx-1)/2)
212 filter_1d%trnsfr(int((filter_1d%nx-1)/2)) = 0.95_rp
213 filter_1d%trnsfr(int((filter_1d%nx-1)/2) + 1) = 0.50_rp
214 filter_1d%trnsfr(int((filter_1d%nx-1)/2) + 2) = 0.05_rp
215 if ((int((filter_1d%nx-1)/2) + 2) .lt. filter_1d%nx) then
216 do i = int((filter_1d%nx-1)/2) + 3, filter_1d%nx
217 filter_1d%trnsfr(i) = 0.0_rp
218 end do
219 end if
220 ! make delta_ratio = (nx-1)/(nt-1) = 2
221 filter_1d%nt = int((filter_1d%nx-1)/2) + 1
222 end if
223
224 call filter_1d%build_1d()
225
226 end subroutine set_ds_filt
227
228end module dynamic_smagorinsky
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Implements the CPU kernel for the smagorinsky_t type.
subroutine, public dynamic_smagorinsky_compute_cpu(if_ext, t, tstep, coef, nut, delta, c_dyn, test_filter, mij, lij, num, den)
Compute eddy viscosity on the CPU.
Implements the device kernel for the smagorinsky_t type.
subroutine, public dynamic_smagorinsky_compute_device(if_ext, t, tstep, coef, nut, delta, c_dyn, test_filter, mij, lij, num, den)
Compute eddy viscosity on the device.
Implements dynamic_smagorinsky_t.
subroutine set_ds_filt(filter_1d)
Set up the test filter.
subroutine dynamic_smagorinsky_free(this)
Destructor for the les_model_t (base) class.
subroutine dynamic_smagorinsky_compute(this, t, tstep)
Compute eddy viscosity.
subroutine dynamic_smagorinsky_init(this, fluid, json)
Constructor.
Implements explicit_filter_t.
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
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:65
integer, parameter, public log_size
Definition log.f90:42
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Utilities.
Definition utils.f90:35
Implements the dynamic Smagorinsky LES model.
Implements the explicit filter for SEM.
Base type of all fluid formulations.
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:65