Neko  0.8.1
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 dirichlet, only : dirichlet_t
40  use utils
41  use, intrinsic :: iso_c_binding, only : c_ptr
42  implicit none
43  private
44 
46  type, public, extends(dirichlet_t) :: facet_normal_t
47  type(coef_t), pointer :: c => null()
48  contains
49  procedure, pass(this) :: apply_scalar => facet_normal_apply_scalar
50  procedure, pass(this) :: apply_vector => facet_normal_apply_vector
51  procedure, pass(this) :: apply_surfvec => facet_normal_apply_surfvec
52  procedure, pass(this) :: apply_surfvec_dev => facet_normal_apply_surfvec_dev
53  end type facet_normal_t
54 
55 contains
56 
58  subroutine facet_normal_apply_scalar(this, x, n, t, tstep)
59  class(facet_normal_t), intent(inout) :: this
60  integer, intent(in) :: n
61  real(kind=rp), intent(inout), dimension(n) :: x
62  real(kind=rp), intent(in), optional :: t
63  integer, intent(in), optional :: tstep
64  end subroutine facet_normal_apply_scalar
65 
67  subroutine facet_normal_apply_vector(this, x, y, z, n, t, tstep)
68  class(facet_normal_t), intent(inout) :: this
69  integer, intent(in) :: n
70  real(kind=rp), intent(inout), dimension(n) :: x
71  real(kind=rp), intent(inout), dimension(n) :: y
72  real(kind=rp), intent(inout), dimension(n) :: z
73  real(kind=rp), intent(in), optional :: t
74  integer, intent(in), optional :: tstep
75  end subroutine facet_normal_apply_vector
76 
78  subroutine facet_normal_apply_surfvec(this, x, y, z, u, v, w, n, t, tstep)
79  class(facet_normal_t), intent(inout) :: this
80  integer, intent(in) :: n
81  real(kind=rp), intent(inout), dimension(n) :: x
82  real(kind=rp), intent(inout), dimension(n) :: y
83  real(kind=rp), intent(inout), dimension(n) :: z
84  real(kind=rp), intent(inout), dimension(n) :: u
85  real(kind=rp), intent(inout), dimension(n) :: v
86  real(kind=rp), intent(inout), dimension(n) :: w
87  real(kind=rp), intent(in), optional :: t
88  integer, intent(in), optional :: tstep
89  integer :: i, m, k, idx(4), facet
90 
91  associate(c => this%coef)
92  m = this%msk(0)
93  do i = 1, m
94  k = this%msk(i)
95  facet = this%facet(i)
96  idx = nonlinear_index(k, c%Xh%lx, c%Xh%lx, c%Xh%lx)
97  select case(facet)
98  case(1,2)
99  x(k) = u(k) * c%nx(idx(2), idx(3), facet, idx(4)) &
100  * c%area(idx(2), idx(3), facet, idx(4))
101  y(k) = v(k) * c%ny(idx(2), idx(3), facet, idx(4)) &
102  * c%area(idx(2), idx(3), facet, idx(4))
103  z(k) = w(k) * c%nz(idx(2), idx(3), facet, idx(4)) &
104  * c%area(idx(2), idx(3), facet, idx(4))
105  case(3,4)
106  x(k) = u(k) * c%nx(idx(1), idx(3), facet, idx(4)) &
107  * c%area(idx(1), idx(3), facet, idx(4))
108  y(k) = v(k) * c%ny(idx(1), idx(3), facet, idx(4)) &
109  * c%area(idx(1), idx(3), facet, idx(4))
110  z(k) = w(k) * c%nz(idx(1), idx(3), facet, idx(4)) &
111  * c%area(idx(1), idx(3), facet, idx(4))
112  case(5,6)
113  x(k) = u(k) * c%nx(idx(1), idx(2), facet, idx(4)) &
114  * c%area(idx(1), idx(2), facet, idx(4))
115  y(k) = v(k) * c%ny(idx(1), idx(2), facet, idx(4)) &
116  * c%area(idx(1), idx(2), facet, idx(4))
117  z(k) = w(k) * c%nz(idx(1), idx(2), facet, idx(4)) &
118  * c%area(idx(1), idx(2), facet, idx(4))
119  end select
120  end do
121  end associate
122 
123  end subroutine facet_normal_apply_surfvec
124 
126  subroutine facet_normal_apply_surfvec_dev(this, x_d, y_d, z_d, &
127  u_d, v_d, w_d, t, tstep)
128  class(facet_normal_t), intent(inout), target :: this
129  type(c_ptr) :: x_d, y_d, z_d, u_d, v_d, w_d
130  real(kind=rp), intent(in), optional :: t
131  integer, intent(in), optional :: tstep
132 
133  associate(c => this%coef)
134  call device_facet_normal_apply_surfvec(this%msk_d, this%facet_d, &
135  x_d, y_d, z_d, u_d, v_d, w_d, &
136  c%nx_d, c%ny_d, c%nz_d, c%area_d, &
137  c%Xh%lx, size(this%msk))
138  end associate
139 
140  end subroutine facet_normal_apply_surfvec_dev
141 
142 end module facet_normal
__device__ void nonlinear_index(const int idx, const int lx, int *index)
Coefficients.
Definition: coef.f90:34
Defines a dirichlet boundary condition.
Definition: dirichlet.f90:34
Dirichlet condition applied in the facet normal direction.
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(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
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
Generic Dirichlet boundary condition on .
Definition: dirichlet.f90:44
Dirichlet condition in facet normal direction.