Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
rough_log_law.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 dofmap, only : dofmap_t
40 use coefs, only : coef_t
42 use wall_model, only : wall_model_t
45 use utils, only : neko_error
46 implicit none
47 private
48
53 type, public, extends(wall_model_t) :: rough_log_law_t
54
56 real(kind=rp) :: kappa = 0.41_rp
58 real(kind=rp) :: b = 0.0_rp
60 real(kind=rp) :: z0 = 0.0_rp
61 contains
63 procedure, pass(this) :: init => rough_log_law_init
65 procedure, pass(this) :: init_from_components => &
68 procedure, pass(this) :: free => rough_log_law_free
70 procedure, pass(this) :: compute => rough_log_law_compute
71 end type rough_log_law_t
72
73contains
81 subroutine rough_log_law_init(this, coef, msk, facet, nu, h_index, json)
82 class(rough_log_law_t), intent(inout) :: this
83 type(coef_t), intent(in) :: coef
84 integer, intent(in) :: msk(:)
85 integer, intent(in) :: facet(:)
86 real(kind=rp), intent(in) :: nu
87 integer, intent(in) :: h_index
88 type(json_file), intent(inout) :: json
89 real(kind=rp) :: kappa, b, z0
90
91 call json_get_or_default(json, "kappa", kappa, 0.41_rp)
92 call json_get(json, "B", b)
93 call json_get(json, "z0", z0)
94
95 call this%init_from_components(coef, msk, facet, nu, h_index, kappa, b, z0)
96 end subroutine rough_log_law_init
97
107 subroutine rough_log_law_init_from_components(this, coef, msk, facet,&
108 nu, h_index, kappa, B, z0)
109 class(rough_log_law_t), intent(inout) :: this
110 type(coef_t), intent(in) :: coef
111 integer, intent(in) :: msk(:)
112 integer, intent(in) :: facet(:)
113 real(kind=rp), intent(in) :: nu
114 integer, intent(in) :: h_index
115 real(kind=rp), intent(in) :: kappa
116 real(kind=rp), intent(in) :: b
117 real(kind=rp), intent(in) :: z0
118
119 if (neko_bcknd_device .eq. 1) then
120 call neko_error("The rough loglaw is only available on the CPU backend.")
121 end if
122
123 call this%init_base(coef, msk, facet, nu, h_index)
124
125 this%kappa = kappa
126 this%B = b
127 this%z0 = z0
128
130
132 subroutine rough_log_law_free(this)
133 class(rough_log_law_t), intent(inout) :: this
134
135 call this%free_base()
136
137 end subroutine rough_log_law_free
138
142 subroutine rough_log_law_compute(this, t, tstep)
143 class(rough_log_law_t), intent(inout) :: this
144 real(kind=rp), intent(in) :: t
145 integer, intent(in) :: tstep
146 type(field_t), pointer :: u
147 type(field_t), pointer :: v
148 type(field_t), pointer :: w
149 integer :: i
150 real(kind=rp) :: ui, vi, wi, magu, utau, normu
151
152 u => neko_field_registry%get_field("u")
153 v => neko_field_registry%get_field("v")
154 w => neko_field_registry%get_field("w")
155
156 do i = 1, this%n_nodes
157 ! Sample the velocity
158 ui = u%x(this%ind_r(i), this%ind_s(i), this%ind_t(i), this%ind_e(i))
159 vi = v%x(this%ind_r(i), this%ind_s(i), this%ind_t(i), this%ind_e(i))
160 wi = w%x(this%ind_r(i), this%ind_s(i), this%ind_t(i), this%ind_e(i))
161
162 ! Project on tangential direction
163 normu = ui * this%n_x%x(i) + vi * this%n_y%x(i) + wi * this%n_z%x(i)
164
165 ui = ui - normu * this%n_x%x(i)
166 vi = vi - normu * this%n_y%x(i)
167 wi = wi - normu * this%n_z%x(i)
168
169 magu = sqrt(ui**2 + vi**2 + wi**2)
170
171 ! Compute the stress
172 utau = (magu - this%B) * this%kappa / log(this%h%x(i) / this%z0)
173
174 ! Distribute according to the velocity vector
175 this%tau_x(i) = -utau**2 * ui / magu
176 this%tau_y(i) = -utau**2 * vi / magu
177 this%tau_z(i) = -utau**2 * wi / magu
178 end do
179
180 end subroutine rough_log_law_compute
181
182
183end module rough_log_law
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
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.
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements rough_log_law_t.
subroutine rough_log_law_init_from_components(this, coef, msk, facet, nu, h_index, kappa, b, z0)
Constructor from components.
subroutine rough_log_law_init(this, coef, msk, facet, nu, h_index, json)
Constructor from JSON.
subroutine rough_log_law_free(this)
Destructor for the rough_log_law_t (base) class.
subroutine rough_log_law_compute(this, t, tstep)
Compute the wall shear stress.
Utilities.
Definition utils.f90:35
Implements wall_model_t.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
Wall model based on the log-law for a rough wall. The formula defining the law is ....
Base abstract type for wall-stress models for wall-modelled LES.