Neko  0.8.99
A portable framework for high-order spectral element flow simulations
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 !
35  use num_types
36  use bc, only : bc_t
37  use, intrinsic :: iso_c_binding, only : c_ptr
38  use utils, only : neko_error, nonlinear_index
39  use coefs, only : coef_t
40  use dirichlet, only : dirichlet_t
41  use neumann, only : neumann_t
42  use math, only : copy
43  use device_math, only : device_copy
44  implicit none
45  private
46 
50  type, public, extends(bc_t) :: shear_stress_t
52  real(kind=rp), allocatable, private :: tau1_(:)
54  real(kind=rp), allocatable, private :: tau2_(:)
58  type(neumann_t) :: neumann1
60  type(neumann_t) :: neumann2
61  contains
62  procedure, pass(this) :: apply_scalar => shear_stress_apply_scalar
63  procedure, pass(this) :: apply_vector => shear_stress_apply_vector
64  procedure, pass(this) :: apply_scalar_dev => shear_stress_apply_scalar_dev
65  procedure, pass(this) :: apply_vector_dev => shear_stress_apply_vector_dev
66  procedure, pass(this) :: init_shear_stress => &
68  procedure, pass(this) :: tau1 => shear_stress_tau1
69  procedure, pass(this) :: tau2 => shear_stress_tau2
70  procedure, pass(this) :: set_stress => shear_stress_set_stress
71  procedure, pass(this) :: shear_stress_finalize
72  procedure, pass(this) :: free => shear_stress_free
73  end type shear_stress_t
74 
75  contains
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%neumann1%apply_scalar(x, n, t, tstep)
104  call this%neumann2%apply_scalar(z, n, t, tstep)
105  call this%dirichlet%apply_scalar(y, 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 
139  subroutine shear_stress_init_shear_stress(this, coef)
140  class(shear_stress_t), intent(inout) :: this
141  type(coef_t), target, intent(in) :: coef
142 
143  call this%init_base(coef)
144 
145  call this%dirichlet%init_base(coef)
146  call this%neumann1%init_base(coef)
147  call this%neumann2%init_base(coef)
148 
149  call this%dirichlet%set_g(0.0_rp)
150 
151  this%coef => coef
152  end subroutine shear_stress_init_shear_stress
153 
156  subroutine shear_stress_finalize(this, tau1, tau2)
157  class(shear_stress_t), intent(inout) :: this
158  real(kind=rp), intent(in) :: tau1
159  real(kind=rp), intent(in) :: tau2
160 
161  call this%neumann1%mark_facets(this%marked_facet)
162  call this%neumann2%mark_facets(this%marked_facet)
163  call this%dirichlet%mark_facets(this%marked_facet)
164 
165  call this%dirichlet%finalize()
166  call this%neumann1%finalize_neumann(tau1)
167  call this%neumann2%finalize_neumann(tau2)
168 
169  end subroutine shear_stress_finalize
170 
172  pure function shear_stress_tau1(this) result(tau1)
173  class(shear_stress_t), intent(in) :: this
174  real(kind=rp) :: tau1(this%msk(0))
175 
176  tau1 = this%tau1_
177  end function shear_stress_tau1
178 
180  pure function shear_stress_tau2(this) result(tau2)
181  class(shear_stress_t), intent(in) :: this
182  real(kind=rp) :: tau2(this%msk(0))
183 
184  tau2 = this%tau2_
185  end function shear_stress_tau2
186 
188  subroutine shear_stress_set_stress(this, tau1, tau2)
189  class(shear_stress_t), intent(inout) :: this
190  real(kind=rp), intent(in) :: tau1(this%msk(0))
191  real(kind=rp), intent(in) :: tau2(this%msk(0))
192 
193  call copy(this%tau1_, tau1, this%msk(0))
194  call copy(this%tau2_, tau2, this%msk(0))
195 
196  end subroutine shear_stress_set_stress
197 
199  subroutine shear_stress_free(this)
200  class(shear_stress_t), target, intent(inout) :: this
201  call this%free_base
202  call this%dirichlet%free
203  call this%neumann1%free
204  call this%neumann2%free
205 
206  end subroutine shear_stress_free
207  end module shear_stress
__device__ void nonlinear_index(const int idx, const int lx, int *index)
Defines a boundary condition.
Definition: bc.f90:34
Coefficients.
Definition: coef.f90:34
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
Defines a dirichlet boundary condition.
Definition: dirichlet.f90:34
Definition: math.f90:60
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:228
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.
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_finalize(this, tau1, tau2)
Finalize by allocating the stress arrays and marking the facets for the bc components.
subroutine shear_stress_free(this)
Destructor.
pure real(kind=rp) function, dimension(this%msk(0)) shear_stress_tau2(this)
Get the stress in the 2nd wall-parallel direction.
subroutine shear_stress_set_stress(this, tau1, tau2)
Set the shear stress components.
subroutine shear_stress_init_shear_stress(this, coef)
Constructor.
subroutine shear_stress_apply_scalar(this, x, n, t, tstep)
Apply shear stress for a scalar field x.
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.
pure real(kind=rp) function, dimension(this%msk(0)) shear_stress_tau1(this)
Get the stress in the 1st wall-parallel direction.
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)
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
Generic Dirichlet boundary condition on .
Definition: dirichlet.f90:44
Generic Neumann boundary condition. This sets the flux of the field to the chosen value.
Definition: neumann.f90:48
A shear stress boundary condition.