Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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, only : rp
37 use math
38 use coefs, only : coef_t
39 use bc, only : bc_t
40 use utils
41 use json_module, only : json_file
42 use, intrinsic :: iso_c_binding, only : c_ptr
43 implicit none
44 private
45
47 type, public, extends(bc_t) :: facet_normal_t
48 contains
49 procedure, pass(this) :: apply_scalar => facet_normal_apply_scalar
50 procedure, pass(this) :: apply_scalar_dev => facet_normal_apply_scalar_dev
51 procedure, pass(this) :: apply_vector => facet_normal_apply_vector
52 procedure, pass(this) :: apply_vector_dev => facet_normal_apply_vector_dev
53 procedure, pass(this) :: apply_surfvec => facet_normal_apply_surfvec
54 procedure, pass(this) :: apply_surfvec_dev => &
57 procedure, pass(this) :: init => facet_normal_init
59 procedure, pass(this) :: init_from_components => &
62 procedure, pass(this) :: free => facet_normal_free
64 procedure, pass(this) :: finalize => facet_normal_finalize
65 end type facet_normal_t
66
67contains
68
72 subroutine facet_normal_init(this, coef, json)
73 class(facet_normal_t), intent(inout), target :: this
74 type(coef_t), intent(in) :: coef
75 type(json_file), intent(inout) ::json
76
77 call this%init_from_components(coef)
78 end subroutine facet_normal_init
79
82 subroutine facet_normal_init_from_components(this, coef)
83 class(facet_normal_t), intent(inout), target :: this
84 type(coef_t), intent(in) :: coef
85
86 call this%init_base(coef)
88
90 subroutine facet_normal_apply_scalar(this, x, n, t, tstep, strong)
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(in), optional :: t
95 integer, intent(in), optional :: tstep
96 logical, intent(in), optional :: strong
97 end subroutine facet_normal_apply_scalar
98
100 subroutine facet_normal_apply_scalar_dev(this, x_d, t, tstep, strong)
101 class(facet_normal_t), intent(inout), target :: this
102 type(c_ptr) :: x_d
103 real(kind=rp), intent(in), optional :: t
104 integer, intent(in), optional :: tstep
105 logical, intent(in), optional :: strong
106
107 end subroutine facet_normal_apply_scalar_dev
108
110 subroutine facet_normal_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, &
111 strong)
112 class(facet_normal_t), intent(inout), target :: this
113 type(c_ptr) :: x_d
114 type(c_ptr) :: y_d
115 type(c_ptr) :: z_d
116 real(kind=rp), intent(in), optional :: t
117 integer, intent(in), optional :: tstep
118 logical, intent(in), optional :: strong
119
120 end subroutine facet_normal_apply_vector_dev
121
123 subroutine facet_normal_apply_vector(this, x, y, z, n, t, tstep, strong)
124 class(facet_normal_t), intent(inout) :: this
125 integer, intent(in) :: n
126 real(kind=rp), intent(inout), dimension(n) :: x
127 real(kind=rp), intent(inout), dimension(n) :: y
128 real(kind=rp), intent(inout), dimension(n) :: z
129 real(kind=rp), intent(in), optional :: t
130 integer, intent(in), optional :: tstep
131 logical, intent(in), optional :: strong
132 end subroutine facet_normal_apply_vector
133
135 subroutine facet_normal_apply_surfvec(this, x, y, z, u, v, w, n, t, tstep)
136 class(facet_normal_t), intent(in) :: this
137 integer, intent(in) :: n
138 real(kind=rp), intent(inout), dimension(n) :: x
139 real(kind=rp), intent(inout), dimension(n) :: y
140 real(kind=rp), intent(inout), dimension(n) :: z
141 real(kind=rp), intent(inout), dimension(n) :: u
142 real(kind=rp), intent(inout), dimension(n) :: v
143 real(kind=rp), intent(inout), dimension(n) :: w
144 real(kind=rp), intent(in), optional :: t
145 integer, intent(in), optional :: tstep
146 integer :: i, m, k, idx(4), facet
147
148 associate(c => this%coef)
149 m = this%msk(0)
150 do i = 1, m
151 k = this%msk(i)
152 facet = this%facet(i)
153 idx = nonlinear_index(k, c%Xh%lx, c%Xh%lx, c%Xh%lx)
154 select case (facet)
155 case (1,2)
156 x(k) = u(k) * c%nx(idx(2), idx(3), facet, idx(4)) &
157 * c%area(idx(2), idx(3), facet, idx(4))
158 y(k) = v(k) * c%ny(idx(2), idx(3), facet, idx(4)) &
159 * c%area(idx(2), idx(3), facet, idx(4))
160 z(k) = w(k) * c%nz(idx(2), idx(3), facet, idx(4)) &
161 * c%area(idx(2), idx(3), facet, idx(4))
162 case (3,4)
163 x(k) = u(k) * c%nx(idx(1), idx(3), facet, idx(4)) &
164 * c%area(idx(1), idx(3), facet, idx(4))
165 y(k) = v(k) * c%ny(idx(1), idx(3), facet, idx(4)) &
166 * c%area(idx(1), idx(3), facet, idx(4))
167 z(k) = w(k) * c%nz(idx(1), idx(3), facet, idx(4)) &
168 * c%area(idx(1), idx(3), facet, idx(4))
169 case (5,6)
170 x(k) = u(k) * c%nx(idx(1), idx(2), facet, idx(4)) &
171 * c%area(idx(1), idx(2), facet, idx(4))
172 y(k) = v(k) * c%ny(idx(1), idx(2), facet, idx(4)) &
173 * c%area(idx(1), idx(2), facet, idx(4))
174 z(k) = w(k) * c%nz(idx(1), idx(2), facet, idx(4)) &
175 * c%area(idx(1), idx(2), facet, idx(4))
176 end select
177 end do
178 end associate
179
180 end subroutine facet_normal_apply_surfvec
181
183 subroutine facet_normal_apply_surfvec_dev(this, x_d, y_d, z_d, &
184 u_d, v_d, w_d, t, tstep)
185 class(facet_normal_t), intent(in), target :: this
186 type(c_ptr) :: x_d, y_d, z_d, u_d, v_d, w_d
187 real(kind=rp), intent(in), optional :: t
188 integer, intent(in), optional :: tstep
189
190 associate(c => this%coef)
191 if (this%msk(0) .gt. 0) then
192 call device_facet_normal_apply_surfvec(this%msk_d, this%facet_d, &
193 x_d, y_d, z_d, u_d, v_d, w_d, c%nx_d, c%ny_d, c%nz_d, c%area_d, &
194 c%Xh%lx, size(this%msk))
195 end if
196 end associate
197
198 end subroutine facet_normal_apply_surfvec_dev
199
201 subroutine facet_normal_free(this)
202 class(facet_normal_t), target, intent(inout) :: this
203
204 call this%free_base()
205
206 end subroutine facet_normal_free
207
209 subroutine facet_normal_finalize(this)
210 class(facet_normal_t), target, intent(inout) :: this
211
212 call this%finalize_base()
213 end subroutine facet_normal_finalize
214
215end 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_init_from_components(this, coef)
Constructor from components.
subroutine facet_normal_apply_scalar_dev(this, x_d, t, tstep, strong)
No-op scalar apply on device.
subroutine facet_normal_init(this, coef, json)
Constructor.
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_finalize(this)
Finalize.
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_vector(this, x, y, z, n, t, tstep, strong)
No-op vector apply.
subroutine facet_normal_apply_scalar(this, x, n, t, tstep, strong)
No-op scalar apply.
subroutine facet_normal_free(this)
Destructor.
subroutine facet_normal_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, strong)
No-op vector apply on device.
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:54
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.