Neko  0.9.99
A portable framework for high-order spectral element flow simulations
facet_normal.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2023, 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
37  use math
38  use coefs, only : coef_t
39  use bc, only : bc_t
40  use utils
41  use, intrinsic :: iso_c_binding, only : c_ptr
42  implicit none
43  private
44 
46  type, public, extends(bc_t) :: facet_normal_t
47  contains
48  procedure, pass(this) :: apply_scalar => facet_normal_apply_scalar
49  procedure, pass(this) :: apply_scalar_dev => facet_normal_apply_scalar_dev
50  procedure, pass(this) :: apply_vector => facet_normal_apply_vector
51  procedure, pass(this) :: apply_vector_dev => facet_normal_apply_vector_dev
52  procedure, pass(this) :: apply_surfvec => facet_normal_apply_surfvec
53  procedure, pass(this) :: apply_surfvec_dev => facet_normal_apply_surfvec_dev
55  procedure, pass(this) :: free => facet_normal_free
56  end type facet_normal_t
57 
58 contains
59 
61  subroutine facet_normal_apply_scalar(this, x, n, t, tstep)
62  class(facet_normal_t), intent(inout) :: this
63  integer, intent(in) :: n
64  real(kind=rp), intent(inout), dimension(n) :: x
65  real(kind=rp), intent(in), optional :: t
66  integer, intent(in), optional :: tstep
67  end subroutine facet_normal_apply_scalar
68 
70  subroutine facet_normal_apply_scalar_dev(this, x_d, t, tstep)
71  class(facet_normal_t), intent(inout), target :: this
72  type(c_ptr) :: x_d
73  real(kind=rp), intent(in), optional :: t
74  integer, intent(in), optional :: tstep
75 
76  end subroutine facet_normal_apply_scalar_dev
77 
79  subroutine facet_normal_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
80  class(facet_normal_t), intent(inout), target :: this
81  type(c_ptr) :: x_d
82  type(c_ptr) :: y_d
83  type(c_ptr) :: z_d
84  real(kind=rp), intent(in), optional :: t
85  integer, intent(in), optional :: tstep
86 
87  end subroutine facet_normal_apply_vector_dev
88 
90  subroutine facet_normal_apply_vector(this, x, y, z, n, t, tstep)
91  class(facet_normal_t), intent(inout) :: this
92  integer, intent(in) :: n
93  real(kind=rp), intent(inout), dimension(n) :: x
94  real(kind=rp), intent(inout), dimension(n) :: y
95  real(kind=rp), intent(inout), dimension(n) :: z
96  real(kind=rp), intent(in), optional :: t
97  integer, intent(in), optional :: tstep
98  end subroutine facet_normal_apply_vector
99 
101  subroutine facet_normal_apply_surfvec(this, x, y, z, u, v, w, n, t, tstep)
102  class(facet_normal_t), intent(inout) :: this
103  integer, intent(in) :: n
104  real(kind=rp), intent(inout), dimension(n) :: x
105  real(kind=rp), intent(inout), dimension(n) :: y
106  real(kind=rp), intent(inout), dimension(n) :: z
107  real(kind=rp), intent(inout), dimension(n) :: u
108  real(kind=rp), intent(inout), dimension(n) :: v
109  real(kind=rp), intent(inout), dimension(n) :: w
110  real(kind=rp), intent(in), optional :: t
111  integer, intent(in), optional :: tstep
112  integer :: i, m, k, idx(4), facet
113 
114  associate(c => this%coef)
115  m = this%msk(0)
116  do i = 1, m
117  k = this%msk(i)
118  facet = this%facet(i)
119  idx = nonlinear_index(k, c%Xh%lx, c%Xh%lx, c%Xh%lx)
120  select case(facet)
121  case(1,2)
122  x(k) = u(k) * c%nx(idx(2), idx(3), facet, idx(4)) &
123  * c%area(idx(2), idx(3), facet, idx(4))
124  y(k) = v(k) * c%ny(idx(2), idx(3), facet, idx(4)) &
125  * c%area(idx(2), idx(3), facet, idx(4))
126  z(k) = w(k) * c%nz(idx(2), idx(3), facet, idx(4)) &
127  * c%area(idx(2), idx(3), facet, idx(4))
128  case(3,4)
129  x(k) = u(k) * c%nx(idx(1), idx(3), facet, idx(4)) &
130  * c%area(idx(1), idx(3), facet, idx(4))
131  y(k) = v(k) * c%ny(idx(1), idx(3), facet, idx(4)) &
132  * c%area(idx(1), idx(3), facet, idx(4))
133  z(k) = w(k) * c%nz(idx(1), idx(3), facet, idx(4)) &
134  * c%area(idx(1), idx(3), facet, idx(4))
135  case(5,6)
136  x(k) = u(k) * c%nx(idx(1), idx(2), facet, idx(4)) &
137  * c%area(idx(1), idx(2), facet, idx(4))
138  y(k) = v(k) * c%ny(idx(1), idx(2), facet, idx(4)) &
139  * c%area(idx(1), idx(2), facet, idx(4))
140  z(k) = w(k) * c%nz(idx(1), idx(2), facet, idx(4)) &
141  * c%area(idx(1), idx(2), facet, idx(4))
142  end select
143  end do
144  end associate
145 
146  end subroutine facet_normal_apply_surfvec
147 
149  subroutine facet_normal_apply_surfvec_dev(this, x_d, y_d, z_d, &
150  u_d, v_d, w_d, t, tstep)
151  class(facet_normal_t), intent(inout), target :: this
152  type(c_ptr) :: x_d, y_d, z_d, u_d, v_d, w_d
153  real(kind=rp), intent(in), optional :: t
154  integer, intent(in), optional :: tstep
155 
156  associate(c => this%coef)
157  call device_facet_normal_apply_surfvec(this%msk_d, this%facet_d, &
158  x_d, y_d, z_d, u_d, v_d, w_d, &
159  c%nx_d, c%ny_d, c%nz_d, c%area_d, &
160  c%Xh%lx, size(this%msk))
161  end associate
162 
163  end subroutine facet_normal_apply_surfvec_dev
164 
166  subroutine facet_normal_free(this)
167  class(facet_normal_t), target, intent(inout) :: this
168 
169  call this%free_base()
170 
171  end subroutine facet_normal_free
172 
173 end module facet_normal
__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
Dirichlet condition applied in the facet normal direction.
subroutine facet_normal_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
No-op vector apply on device.
subroutine facet_normal_apply_surfvec_dev(this, x_d, y_d, z_d, u_d, v_d, w_d, t, tstep)
Apply in facet normal direction (vector valued, device version)
subroutine facet_normal_apply_surfvec(this, x, y, z, u, v, w, n, t, tstep)
Apply in facet normal direction (vector valued)
subroutine facet_normal_apply_scalar_dev(this, x_d, t, tstep)
No-op scalar apply on device.
subroutine facet_normal_free(this)
Destructor.
subroutine facet_normal_apply_scalar(this, x, n, t, tstep)
No-op scalar apply.
subroutine facet_normal_apply_vector(this, x, y, z, n, t, tstep)
No-op vector apply.
Definition: math.f90:60
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
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
Dirichlet condition in facet normal direction.