Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
shear_stress.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!
36 use num_types, only : rp
37 use bc, only : bc_t
38 use, intrinsic :: iso_c_binding, only : c_ptr
39 use utils, only : neko_error
40 use coefs, only : coef_t
41 use symmetry, only : symmetry_t
42 use neumann, only : neumann_t
43 implicit none
44 private
45
48 type, public, extends(bc_t) :: shear_stress_t
49 ! This bc takes care of setting the wall-normal component to zero.
50 ! It can be passed to associated bc lists, which take care of masking
51 ! changes in residuals and solution increments.
53
55 type(neumann_t) :: neumann_x
57 type(neumann_t) :: neumann_y
59 type(neumann_t) :: neumann_z
60 contains
61 procedure, pass(this) :: apply_scalar => shear_stress_apply_scalar
62 procedure, pass(this) :: apply_vector => shear_stress_apply_vector
63 procedure, pass(this) :: apply_scalar_dev => shear_stress_apply_scalar_dev
64 procedure, pass(this) :: apply_vector_dev => shear_stress_apply_vector_dev
65 procedure, pass(this) :: init_shear_stress => &
67 procedure, pass(this) :: set_stress_scalar => &
69 procedure, pass(this) :: set_stress_array => &
71 generic :: set_stress => set_stress_scalar, set_stress_array
72 procedure, pass(this) :: free => shear_stress_free
73 end type shear_stress_t
74
75contains
76
78 subroutine shear_stress_apply_scalar(this, x, n, t, tstep)
79 class(shear_stress_t), intent(inout) :: this
80 integer, intent(in) :: n
81 real(kind=rp), intent(inout), dimension(n) :: x
82 real(kind=rp), intent(in), optional :: t
83 integer, intent(in), optional :: tstep
84 integer :: i, m, k, facet
85 ! Store non-linear index
86 integer :: idx(4)
87
88 call neko_error("The shear stress bc is not applicable to scalar fields.")
89
90 end subroutine shear_stress_apply_scalar
91
94 subroutine shear_stress_apply_vector(this, x, y, z, n, t, tstep)
95 class(shear_stress_t), intent(inout) :: this
96 integer, intent(in) :: n
97 real(kind=rp), intent(inout), dimension(n) :: x
98 real(kind=rp), intent(inout), dimension(n) :: y
99 real(kind=rp), intent(inout), dimension(n) :: z
100 real(kind=rp), intent(in), optional :: t
101 integer, intent(in), optional :: tstep
102
103 call this%neumann_x%apply_scalar(x, n, t, tstep)
104 call this%neumann_y%apply_scalar(y, n, t, tstep)
105 call this%neumann_z%apply_scalar(z, n, t, tstep)
106
107 end subroutine shear_stress_apply_vector
108
111 subroutine shear_stress_apply_scalar_dev(this, x_d, t, tstep)
112 class(shear_stress_t), intent(inout), target :: this
113 type(c_ptr) :: x_d
114 real(kind=rp), intent(in), optional :: t
115 integer, intent(in), optional :: tstep
116
117 call neko_error("shear_stress bc not implemented on the device")
118
119 end subroutine shear_stress_apply_scalar_dev
120
123 subroutine shear_stress_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
124 class(shear_stress_t), intent(inout), target :: this
125 type(c_ptr) :: x_d
126 type(c_ptr) :: y_d
127 type(c_ptr) :: z_d
128 real(kind=rp), intent(in), optional :: t
129 integer, intent(in), optional :: tstep
130
131 call neko_error("shear_stress bc not implemented on the device")
132
133 end subroutine shear_stress_apply_vector_dev
134
138 subroutine shear_stress_init_shear_stress(this, coef)
139 class(shear_stress_t), intent(inout) :: this
140 type(coef_t), target, intent(in) :: coef
141
142 this%coef => coef
143
144 call this%symmetry%free()
145 call this%symmetry%init_base(coef)
146 call this%symmetry%mark_facets(this%marked_facet)
147 call this%symmetry%finalize()
148 call this%symmetry%init(coef)
149
150 call this%neumann_x%init_base(coef)
151 call this%neumann_y%init_base(coef)
152 call this%neumann_z%init_base(coef)
153
154 call this%neumann_x%mark_facets(this%marked_facet)
155 call this%neumann_y%mark_facets(this%marked_facet)
156 call this%neumann_z%mark_facets(this%marked_facet)
157
158 call this%neumann_x%finalize_neumann(0.0_rp)
159 call this%neumann_y%finalize_neumann(0.0_rp)
160 call this%neumann_z%finalize_neumann(0.0_rp)
161
162
163 end subroutine shear_stress_init_shear_stress
164
166 subroutine shear_stress_set_stress_scalar(this, tau_x, tau_y, tau_z)
167 class(shear_stress_t), intent(inout) :: this
168 real(kind=rp), intent(in) :: tau_x
169 real(kind=rp), intent(in) :: tau_y
170 real(kind=rp), intent(in) :: tau_z
171
172 ! Calls finalize and allocates the flux arrays
173 call this%neumann_x%set_flux(tau_x)
174 call this%neumann_y%set_flux(tau_y)
175 call this%neumann_z%set_flux(tau_z)
176
177
178 end subroutine shear_stress_set_stress_scalar
179
181 subroutine shear_stress_set_stress_array(this, tau_x, tau_y, tau_z)
182 class(shear_stress_t), intent(inout) :: this
183 real(kind=rp), intent(in) :: tau_x(this%msk(0))
184 real(kind=rp), intent(in) :: tau_y(this%msk(0))
185 real(kind=rp), intent(in) :: tau_z(this%msk(0))
186
187 call this%neumann_x%set_flux(tau_x)
188 call this%neumann_y%set_flux(tau_y)
189 call this%neumann_z%set_flux(tau_z)
190
191 end subroutine shear_stress_set_stress_array
192
194 subroutine shear_stress_free(this)
195 class(shear_stress_t), target, intent(inout) :: this
196 call this%free_base
197 call this%symmetry%free
198
199 call this%neumann_x%free
200 call this%neumann_y%free
201 call this%neumann_z%free
202
203 end subroutine shear_stress_free
204end module shear_stress
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
Defines a Neumann boundary condition.
Definition neumann.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a shear stress boundary condition for a vector field. Maintainer: Timofey Mukha.
subroutine shear_stress_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
Boundary condition apply for a generic shear_stress condition to vectors x, y and z (device version)
subroutine shear_stress_free(this)
Destructor.
subroutine shear_stress_init_shear_stress(this, coef)
Additional constructor that should be run after the finalization of the bc. Similar to the symmetry c...
subroutine shear_stress_set_stress_array(this, tau_x, tau_y, tau_z)
Set the shear stress components.
subroutine shear_stress_apply_scalar(this, x, n, t, tstep)
Apply shear stress for a scalar field x.
subroutine shear_stress_set_stress_scalar(this, tau_x, tau_y, tau_z)
Set the value of the shear stress vector using 3 scalars.
subroutine shear_stress_apply_vector(this, x, y, z, n, t, tstep)
Boundary condition apply for a generic shear_stress condition to vectors x, y and z.
subroutine shear_stress_apply_scalar_dev(this, x_d, t, tstep)
Boundary condition apply for a generic shear_stress condition to a vector x (device version)
Mixed Dirichlet-Neumann axis aligned symmetry plane.
Definition symmetry.f90:34
Utilities.
Definition utils.f90:35
Base type for a boundary condition.
Definition bc.f90:51
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
A Neumann boundary condition for scalar fields. Sets the flux of the field to the chosen value.
Definition neumann.f90:48
A shear stress boundary condition.
Mixed Dirichlet-Neumann symmetry plane condition.
Definition symmetry.f90:46