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 contains
67 procedure, pass(this) :: init => spalding_init
70 procedure, pass(this) :: partial_init => spalding_partial_init
72 procedure, pass(this) :: finalize => spalding_finalize
74 procedure, pass(this) :: init_from_components => &
77 procedure, pass(this) :: free => spalding_free
79 procedure, pass(this) :: compute_nu => spalding_compute_nu
81 procedure, pass(this) :: compute => spalding_compute
82 end type spalding_t
83
84contains
92 subroutine spalding_init(this, scheme_name, coef, msk, facet, h_index, json)
93 class(spalding_t), intent(inout) :: this
94 character(len=*), intent(in) :: scheme_name
95 type(coef_t), intent(in) :: coef
96 integer, intent(in) :: msk(:)
97 integer, intent(in) :: facet(:)
98 integer, intent(in) :: h_index
99 type(json_file), intent(inout) :: json
100 real(kind=rp) :: kappa, b
101
102 call json_get_or_lookup(json, "kappa", kappa)
103 call json_get_or_lookup(json, "B", b)
104
105 call this%init_from_components(scheme_name, coef, msk, facet, h_index, &
106 kappa, b)
107 end subroutine spalding_init
108
112 subroutine spalding_partial_init(this, coef, json)
113 class(spalding_t), intent(inout) :: this
114 type(coef_t), intent(in) :: coef
115 type(json_file), intent(inout) :: json
116 character(len=LOG_SIZE) :: log_buf
117
118 call this%partial_init_base(coef, json)
119 call json_get_or_lookup(json, "kappa", this%kappa)
120 call json_get_or_lookup(json, "B", this%B)
121
122 call neko_log%section('Wall model')
123 write(log_buf, '(A)') 'Model : Spalding'
124 call neko_log%message(log_buf)
125 write(log_buf, '(A, E15.7)') 'kappa : ', this%kappa
126 call neko_log%message(log_buf)
127 write(log_buf, '(A, E15.7)') 'B : ', this%B
128 call neko_log%message(log_buf)
129 call neko_log%end_section()
130
131 end subroutine spalding_partial_init
132
136 subroutine spalding_finalize(this, msk, facet)
137 class(spalding_t), intent(inout) :: this
138 integer, intent(in) :: msk(:)
139 integer, intent(in) :: facet(:)
140
141 call this%finalize_base(msk, facet)
142 call this%nu%init(this%n_nodes)
143 end subroutine spalding_finalize
144
153 subroutine spalding_init_from_components(this, scheme_name, coef, msk, &
154 facet, h_index, kappa, B)
155 class(spalding_t), intent(inout) :: this
156 character(len=*), intent(in) :: scheme_name
157 type(coef_t), intent(in) :: coef
158 integer, intent(in) :: msk(:)
159 integer, intent(in) :: facet(:)
160 integer, intent(in) :: h_index
161 real(kind=rp), intent(in) :: kappa
162 real(kind=rp), intent(in) :: b
163
164 call this%free()
165 call this%init_base(scheme_name, coef, msk, facet, h_index)
166
167 this%kappa = kappa
168 this%B = b
169
170 call this%nu%init(this%n_nodes)
171 end subroutine spalding_init_from_components
172
174 subroutine spalding_compute_nu(this)
175 class(spalding_t), intent(inout) :: this
176 type(field_t), pointer :: temp
177 integer :: idx
178
179 call neko_scratch_registry%request_field(temp, idx, .false.)
180 call field_invcol3(temp, this%mu, this%rho)
181
182 if (neko_bcknd_device .eq. 1) then
183 call device_masked_gather_copy_0(this%nu%x_d, temp%x_d, this%msk_d, &
184 temp%size(), this%nu%size())
185 else
186 call masked_gather_copy_0(this%nu%x, temp%x, this%msk, temp%size(), &
187 this%nu%size())
188 end if
189
190 call neko_scratch_registry%relinquish_field(idx)
191 end subroutine spalding_compute_nu
192
194 subroutine spalding_free(this)
195 class(spalding_t), intent(inout) :: this
196
197 call this%free_base()
198
199 end subroutine spalding_free
200
204 subroutine spalding_compute(this, t, tstep)
205 class(spalding_t), intent(inout) :: this
206 real(kind=rp), intent(in) :: t
207 integer, intent(in) :: tstep
208 type(field_t), pointer :: u
209 type(field_t), pointer :: v
210 type(field_t), pointer :: w
211 integer :: i
212 real(kind=rp) :: ui, vi, wi, magu, utau, normu, guess
213
214 call this%compute_nu()
215
216 u => neko_registry%get_field("u")
217 v => neko_registry%get_field("v")
218 w => neko_registry%get_field("w")
219
220 if (neko_bcknd_device .eq. 1) then
221 call spalding_compute_device(u%x_d, v%x_d, w%x_d, this%ind_r_d, &
222 this%ind_s_d, this%ind_t_d, this%ind_e_d, &
223 this%n_x%x_d, this%n_y%x_d, this%n_z%x_d, &
224 this%nu%x_d, this%h%x_d, &
225 this%tau_x%x_d, this%tau_y%x_d, this%tau_z%x_d, &
226 this%n_nodes, u%Xh%lx, &
227 this%kappa, this%B, tstep)
228 else
229 call spalding_compute_cpu(u%x, v%x, w%x, &
230 this%ind_r, this%ind_s, this%ind_t, this%ind_e, &
231 this%n_x%x, this%n_y%x, this%n_z%x, &
232 this%nu%x, this%h%x, &
233 this%tau_x%x, this%tau_y%x, this%tau_z%x, &
234 this%n_nodes, u%Xh%lx, u%msh%nelv, &
235 this%kappa, this%B, tstep)
236 end if
237
238 end subroutine spalding_compute
239end 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__ 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:313
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, 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, 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:155
subroutine spalding_partial_init(this, coef, json)
Constructor from JSON.
Definition spalding.f90:113
subroutine spalding_compute_nu(this)
Compute the kinematic viscosity vector.
Definition spalding.f90:175
subroutine spalding_free(this)
Destructor for the spalding_t (base) class.
Definition spalding.f90:195
subroutine spalding_finalize(this, msk, facet)
Finalize the construction using the mask and facet arrays of the bc.
Definition spalding.f90:137
subroutine spalding_init(this, scheme_name, coef, msk, facet, h_index, json)
Constructor from JSON.
Definition spalding.f90:93
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:62
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.