Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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 math
37 use field_list, only : field_list_t
38 use field, only : field_t
39 use les_model, only : les_model_t
40 use dofmap , only : dofmap_t
42 use json_module, only : json_file
43 use utils, only : neko_error
45 use coefs, only : coef_t
48 use logger, only : log_size, neko_log
50 implicit none
51 private
52
55 type, public, extends(les_model_t) :: dynamic_smagorinsky_t
57 type(field_t) :: c_dyn
59 type(elementwise_filter_t) :: test_filter
63 type(field_t) :: mij(6)
65 type(field_t) :: lij(6)
67 type(field_t) :: num
69 type(field_t) :: den
70 contains
72 procedure, pass(this) :: init => dynamic_smagorinsky_init
74 procedure, pass(this) :: free => dynamic_smagorinsky_free
76 procedure, pass(this) :: compute => dynamic_smagorinsky_compute
77
79
80contains
85 subroutine dynamic_smagorinsky_init(this, dofmap, coef, json)
86 class(dynamic_smagorinsky_t), intent(inout) :: this
87 type(dofmap_t), intent(in) :: dofmap
88 type(coef_t), intent(in) :: coef
89 type(json_file), intent(inout) :: json
90 character(len=:), allocatable :: nut_name !! The name of the SGS viscosity field.
91 integer :: i
92 character(len=:), allocatable :: delta_type
93 character(len=:), allocatable :: test_filter_type
94 character(len=LOG_SIZE) :: log_buf
95
96 call json_get_or_default(json, "nut_field", nut_name, "nut")
97 call json_get_or_default(json, "delta_type", delta_type, "pointwise")
98 call json_get_or_default(json, "test_filter_type", test_filter_type, "nonBoyd")
99
100 call this%free()
101 call this%init_base(dofmap, coef, nut_name, delta_type)
102 ! Filter assumes lx = ly = lz
103 call this%test_filter%init(dofmap%xh%lx, test_filter_type)
104 call set_ds_filt(this%test_filter)
105
106 call neko_log%section('LES model')
107 write(log_buf, '(A)') 'Model : Dynamic Smagorinsky'
108 call neko_log%message(log_buf)
109 write(log_buf, '(A, A)') 'Delta evaluation : ', delta_type
110 call neko_log%message(log_buf)
111 write(log_buf, '(A, A)') 'Test filter type : ', test_filter_type
112 call neko_log%message(log_buf)
113 call neko_log%end_section()
114
115 call this%c_dyn%init(dofmap, "ds_c_dyn")
116 call this%num%init(dofmap, "ds_num")
117 call this%den%init(dofmap, "ds_den")
118
119 do i = 1, 6
120 call this%mij(i)%init(dofmap)
121 call this%lij(i)%init(dofmap)
122 end do
123
124 end subroutine dynamic_smagorinsky_init
125
128 class(dynamic_smagorinsky_t), intent(inout) :: this
129 integer :: i
130
131 call this%c_dyn%free()
132 do i = 1, 6
133 call this%mij(i)%free()
134 call this%lij(i)%free()
135 end do
136 call this%num%free()
137 call this%den%free()
138 call this%test_filter%free()
139 call this%free_base()
140
141 end subroutine dynamic_smagorinsky_free
142
146 subroutine dynamic_smagorinsky_compute(this, t, tstep)
147 class(dynamic_smagorinsky_t), intent(inout) :: this
148 real(kind=rp), intent(in) :: t
149 integer, intent(in) :: tstep
150
151 if (neko_bcknd_device .eq. 1) then
152 call dynamic_smagorinsky_compute_device(t, tstep, this%coef, this%nut, &
153 this%delta, this%c_dyn, this%test_filter, &
154 this%mij, this%lij, this%num, this%den)
155 else
156 call dynamic_smagorinsky_compute_cpu(t, tstep, this%coef, this%nut, &
157 this%delta, this%c_dyn, this%test_filter, &
158 this%mij, this%lij, this%num, this%den)
159 end if
160
161 end subroutine dynamic_smagorinsky_compute
162
164 subroutine set_ds_filt(filter_1d)
165 type(elementwise_filter_t), intent(inout) :: filter_1d
166 integer :: i
167
168 if (filter_1d%nx .le. 2) then
169 call neko_error("Dynamic Smagorinsky model error: test filter is not &
170 &defined for the current polynomial order")
171 end if
172 if (mod(filter_1d%nx,2) .eq. 0) then ! number of grid spacing is odd
173 ! cutoff at polynomial order int((filter_1d%nx)/2)
174 filter_1d%trnsfr(int((filter_1d%nx)/2)) = 0.95_rp
175 filter_1d%trnsfr(int((filter_1d%nx)/2) + 1) = 0.50_rp
176 filter_1d%trnsfr(int((filter_1d%nx)/2) + 2) = 0.05_rp
177 if ((int((filter_1d%nx)/2) + 2) .lt. filter_1d%nx) then
178 do i = int((filter_1d%nx)/2) + 3, filter_1d%nx
179 filter_1d%trnsfr(i) = 0.0_rp
180 end do
181 end if
182 ! make delta_ratio = (nx-1)/(nt-1) as close to 2
183 filter_1d%nt = int(filter_1d%nx/2) + 1
184 else ! number of grid spacing is even
185 ! cutoff at polynomial order int((filter_1d%nx-1)/2)
186 filter_1d%trnsfr(int((filter_1d%nx-1)/2)) = 0.95_rp
187 filter_1d%trnsfr(int((filter_1d%nx-1)/2) + 1) = 0.50_rp
188 filter_1d%trnsfr(int((filter_1d%nx-1)/2) + 2) = 0.05_rp
189 if ((int((filter_1d%nx-1)/2) + 2) .lt. filter_1d%nx) then
190 do i = int((filter_1d%nx-1)/2) + 3, filter_1d%nx
191 filter_1d%trnsfr(i) = 0.0_rp
192 end do
193 end if
194 ! make delta_ratio = (nx-1)/(nt-1) = 2
195 filter_1d%nt = int((filter_1d%nx-1)/2) + 1
196 end if
197
198 call filter_1d%build_1d()
199
200 end subroutine set_ds_filt
201
202end module dynamic_smagorinsky
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.
Coefficients.
Definition coef.f90:34
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(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(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_init(this, dofmap, coef, json)
Constructor.
subroutine dynamic_smagorinsky_free(this)
Destructor for the les_model_t (base) class.
subroutine dynamic_smagorinsky_compute(this, t, tstep)
Compute eddy viscosity.
Implements explicit_filter_t.
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
Definition math.f90:60
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
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Implements the dynamic Smagorinsky LES model.
Implements the explicit filter for SEM.
field_list_t, To be able to group fields together
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:51