Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
spalding.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!
33!
36 use field, only: field_t
37 use num_types, only : rp
38 use json_module, only : json_file
39 use coefs, only : coef_t
41 use wall_model, only : wall_model_t
42 use registry, only : neko_registry
46 use field_math, only: field_invcol3
47 use vector, only : vector_t
48 use math, only: masked_gather_copy_0
51 use logger, only : log_size, neko_log
52
53 implicit none
54 private
55
58 type, public, extends(wall_model_t) :: spalding_t
60 real(kind=rp) :: kappa = 0.41_rp
62 real(kind=rp) :: b = 5.2_rp
64 type(vector_t) :: nu
65 ! The fluid density at the boundary
66 type(vector_t) :: rho_w
67 contains
69 procedure, pass(this) :: init => spalding_init
72 procedure, pass(this) :: partial_init => spalding_partial_init
74 procedure, pass(this) :: finalize => spalding_finalize
76 procedure, pass(this) :: init_from_components => &
79 procedure, pass(this) :: free => spalding_free
81 procedure, pass(this) :: compute_nu => spalding_compute_nu
83 procedure, pass(this) :: compute => spalding_compute
84 end type spalding_t
85
86contains
94 subroutine spalding_init(this, scheme_name, coef, msk, facet, h_index, json)
95 class(spalding_t), intent(inout) :: this
96 character(len=*), intent(in) :: scheme_name
97 type(coef_t), intent(in) :: coef
98 integer, intent(in) :: msk(:)
99 integer, intent(in) :: facet(:)
100 integer, intent(in) :: h_index
101 type(json_file), intent(inout) :: json
102 real(kind=rp) :: kappa, b
103
104 call json_get_or_lookup(json, "kappa", kappa)
105 call json_get_or_lookup(json, "B", b)
106
107 call this%init_from_components(scheme_name, coef, msk, facet, h_index, &
108 kappa, b)
109 end subroutine spalding_init
110
114 subroutine spalding_partial_init(this, coef, json)
115 class(spalding_t), intent(inout) :: this
116 type(coef_t), intent(in) :: coef
117 type(json_file), intent(inout) :: json
118 character(len=LOG_SIZE) :: log_buf
119
120 call this%partial_init_base(coef, json)
121 call json_get_or_lookup(json, "kappa", this%kappa)
122 call json_get_or_lookup(json, "B", this%B)
123
124 call neko_log%section('Wall model')
125 write(log_buf, '(A)') 'Model : Spalding'
126 call neko_log%message(log_buf)
127 write(log_buf, '(A, E15.7)') 'kappa : ', this%kappa
128 call neko_log%message(log_buf)
129 write(log_buf, '(A, E15.7)') 'B : ', this%B
130 call neko_log%message(log_buf)
131 call neko_log%end_section()
132
133 end subroutine spalding_partial_init
134
138 subroutine spalding_finalize(this, msk, facet)
139 class(spalding_t), intent(inout) :: this
140 integer, intent(in) :: msk(:)
141 integer, intent(in) :: facet(:)
142
143 call this%finalize_base(msk, facet)
144 call this%nu%init(this%n_nodes)
145 call this%rho_w%init(this%n_nodes)
146 end subroutine spalding_finalize
147
156 subroutine spalding_init_from_components(this, scheme_name, coef, msk, &
157 facet, h_index, kappa, B)
158 class(spalding_t), intent(inout) :: this
159 character(len=*), intent(in) :: scheme_name
160 type(coef_t), intent(in) :: coef
161 integer, intent(in) :: msk(:)
162 integer, intent(in) :: facet(:)
163 integer, intent(in) :: h_index
164 real(kind=rp), intent(in) :: kappa
165 real(kind=rp), intent(in) :: b
166
167 call this%free()
168 call this%init_base(scheme_name, coef, msk, facet, h_index)
169
170 this%kappa = kappa
171 this%B = b
172
173 call this%nu%init(this%n_nodes)
174 call this%rho_w%init(this%n_nodes)
175 end subroutine spalding_init_from_components
176
178 subroutine spalding_compute_nu(this)
179 class(spalding_t), intent(inout) :: this
180 type(field_t), pointer :: temp
181 integer :: idx
182
183 call neko_scratch_registry%request_field(temp, idx, .false.)
184 call field_invcol3(temp, this%mu, this%rho)
185
186 if (neko_bcknd_device .eq. 1) then
187 call device_masked_gather_copy_0(this%nu%x_d, temp%x_d, this%msk_d, &
188 temp%size(), this%nu%size())
189 call device_masked_gather_copy_0(this%rho_w%x_d, this%rho%x_d, this%msk_d, &
190 this%rho%size(), this%rho_w%size())
191 else
192 call masked_gather_copy_0(this%nu%x, temp%x, this%msk, temp%size(), &
193 this%nu%size())
194 call masked_gather_copy_0(this%rho_w%x, this%rho%x, this%msk, &
195 this%rho%size(), this%rho_w%size())
196 end if
197
198 call neko_scratch_registry%relinquish_field(idx)
199 end subroutine spalding_compute_nu
200
202 subroutine spalding_free(this)
203 class(spalding_t), intent(inout) :: this
204
205 call this%rho_w%free()
206 call this%free_base()
207
208 end subroutine spalding_free
209
213 subroutine spalding_compute(this, t, tstep)
214 class(spalding_t), intent(inout) :: this
215 real(kind=rp), intent(in) :: t
216 integer, intent(in) :: tstep
217 type(field_t), pointer :: u
218 type(field_t), pointer :: v
219 type(field_t), pointer :: w
220 integer :: i
221 real(kind=rp) :: ui, vi, wi, magu, utau, normu, guess
222
223 call this%compute_nu()
224
225 u => neko_registry%get_field("u")
226 v => neko_registry%get_field("v")
227 w => neko_registry%get_field("w")
228
229 if (neko_bcknd_device .eq. 1) then
230 call spalding_compute_device(u%x_d, v%x_d, w%x_d, this%ind_r_d, &
231 this%ind_s_d, this%ind_t_d, this%ind_e_d, &
232 this%n_x%x_d, this%n_y%x_d, this%n_z%x_d, &
233 this%nu%x_d, this%rho_w%x_d, this%h%x_d, &
234 this%tau_x%x_d, this%tau_y%x_d, this%tau_z%x_d, &
235 this%n_nodes, u%Xh%lx, &
236 this%kappa, this%B, tstep)
237 else
238 call spalding_compute_cpu(u%x, v%x, w%x, &
239 this%ind_r, this%ind_s, this%ind_t, this%ind_e, &
240 this%n_x%x, this%n_y%x, this%n_z%x, &
241 this%nu%x, this%rho_w%x, this%h%x, &
242 this%tau_x%x, this%tau_y%x, this%tau_z%x, &
243 this%n_nodes, u%Xh%lx, u%msh%nelv, &
244 this%kappa, this%B, tstep)
245 end if
246
247 end subroutine spalding_compute
248end module spalding
__global__ void spalding_compute(const T *__restrict__ u_d, const T *__restrict__ v_d, const T *__restrict__ w_d, const int *__restrict__ ind_r_d, const int *__restrict__ ind_s_d, const int *__restrict__ ind_t_d, const int *__restrict__ ind_e_d, const T *__restrict__ n_x_d, const T *__restrict__ n_y_d, const T *__restrict__ n_z_d, const T *__restrict__ nu_d, const T *__restrict__ rho_w_d, const T *__restrict__ h_d, T *__restrict__ tau_x_d, T *__restrict__ tau_y_d, T *__restrict__ tau_z_d, const int n_nodes, const int lx, const T kappa, const T B, const int tstep)
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Coefficients.
Definition coef.f90:34
subroutine, public device_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
subroutine, public field_invcol3(a, b, c, n)
Invert a vector .
Defines a field.
Definition field.f90:34
Utilities for retrieving parameters from the case files.
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
Definition math.f90:60
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:315
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
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Implements the CPU kernel for the spalding_t type.
subroutine, public spalding_compute_cpu(u, v, w, ind_r, ind_s, ind_t, ind_e, n_x, n_y, n_z, nu, rho_w, h, tau_x, tau_y, tau_z, n_nodes, lx, nelv, kappa, b, tstep)
Compute the wall shear stress on cpu using Spalding's model.
Implements the device kernel for the spalding_t type.
subroutine, public spalding_compute_device(u_d, v_d, w_d, ind_r_d, ind_s_d, ind_t_d, ind_e_d, n_x_d, n_y_d, n_z_d, nu_d, rho_w_d, h_d, tau_x_d, tau_y_d, tau_z_d, n_nodes, lx, kappa, b, tstep)
Compute the wall shear stress on device using Spalding's model.
Implements spalding_t.
Definition spalding.f90:35
subroutine spalding_init_from_components(this, scheme_name, coef, msk, facet, h_index, kappa, b)
Constructor from components.
Definition spalding.f90:158
subroutine spalding_partial_init(this, coef, json)
Constructor from JSON.
Definition spalding.f90:115
subroutine spalding_compute_nu(this)
Compute the kinematic viscosity vector.
Definition spalding.f90:179
subroutine spalding_free(this)
Destructor for the spalding_t (base) class.
Definition spalding.f90:203
subroutine spalding_finalize(this, msk, facet)
Finalize the construction using the mask and facet arrays of the bc.
Definition spalding.f90:139
subroutine spalding_init(this, scheme_name, coef, msk, facet, h_index, json)
Constructor from JSON.
Definition spalding.f90:95
Defines a vector.
Definition vector.f90:34
Implements wall_model_t.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Wall model based on Spalding's law of the wall. Reference: http://dx.doi.org/10.1115/1....
Definition spalding.f90:58
Base abstract type for wall-stress models for wall-modelled LES.