Neko 1.99.2
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-2025, 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, only : rp
37 use math, only : cfill_mask
40 use vector, only : vector_t
41 use coefs, only : coef_t
42 use bc, only : bc_t
44 use json_module, only : json_file
45 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
46 use htable, only : htable_i4_t
49 use time_state, only : time_state_t
50 implicit none
51 private
52
54 type, public, extends(bc_t) :: facet_normal_t
55 integer, allocatable :: unique_mask(:)
56 type(c_ptr) :: unique_mask_d = c_null_ptr
57 type(vector_t) :: nx, ny, nz, work
58 contains
59 procedure, pass(this) :: apply_scalar => facet_normal_apply_scalar
60 procedure, pass(this) :: apply_scalar_dev => facet_normal_apply_scalar_dev
61 procedure, pass(this) :: apply_vector => facet_normal_apply_vector
62 procedure, pass(this) :: apply_vector_dev => facet_normal_apply_vector_dev
63 procedure, pass(this) :: apply_surfvec => facet_normal_apply_surfvec
64 procedure, pass(this) :: apply_surfvec_dev => &
67 procedure, pass(this) :: init => facet_normal_init
69 procedure, pass(this) :: init_from_components => &
72 procedure, pass(this) :: free => facet_normal_free
74 procedure, pass(this) :: finalize => facet_normal_finalize
75 end type facet_normal_t
76
77contains
78
82 subroutine facet_normal_init(this, coef, json)
83 class(facet_normal_t), intent(inout), target :: this
84 type(coef_t), target, intent(in) :: coef
85 type(json_file), intent(inout) :: json
86
87 call this%init_from_components(coef)
88 end subroutine facet_normal_init
89
92 subroutine facet_normal_init_from_components(this, coef)
93 class(facet_normal_t), intent(inout), target :: this
94 type(coef_t), target, intent(in) :: coef
95
96 call this%init_base(coef)
98
100 subroutine facet_normal_apply_scalar(this, x, n, time, strong)
101 class(facet_normal_t), intent(inout) :: this
102 integer, intent(in) :: n
103 real(kind=rp), intent(inout), dimension(n) :: x
104 type(time_state_t), intent(in), optional :: time
105 logical, intent(in), optional :: strong
106 end subroutine facet_normal_apply_scalar
107
109 subroutine facet_normal_apply_scalar_dev(this, x_d, time, strong, strm)
110 class(facet_normal_t), intent(inout), target :: this
111 type(c_ptr), intent(inout) :: x_d
112 type(time_state_t), intent(in), optional :: time
113 logical, intent(in), optional :: strong
114 type(c_ptr), intent(inout) :: strm
115
116 end subroutine facet_normal_apply_scalar_dev
117
119 subroutine facet_normal_apply_vector_dev(this, x_d, y_d, z_d, time, &
120 strong, strm)
121 class(facet_normal_t), intent(inout), target :: this
122 type(c_ptr), intent(inout) :: x_d
123 type(c_ptr), intent(inout) :: y_d
124 type(c_ptr), intent(inout) :: z_d
125 type(time_state_t), intent(in), optional :: time
126 logical, intent(in), optional :: strong
127 type(c_ptr), intent(inout) :: strm
128
129 end subroutine facet_normal_apply_vector_dev
130
132 subroutine facet_normal_apply_vector(this, x, y, z, n, time, strong)
133 class(facet_normal_t), intent(inout) :: this
134 integer, intent(in) :: n
135 real(kind=rp), intent(inout), dimension(n) :: x
136 real(kind=rp), intent(inout), dimension(n) :: y
137 real(kind=rp), intent(inout), dimension(n) :: z
138 type(time_state_t), intent(in), optional :: time
139 logical, intent(in), optional :: strong
140 end subroutine facet_normal_apply_vector
141
143 subroutine facet_normal_apply_surfvec(this, x, y, z, u, v, w, n, time)
144 class(facet_normal_t), intent(in) :: this
145 integer, intent(in) :: n
146 real(kind=rp), intent(inout), dimension(n) :: x
147 real(kind=rp), intent(inout), dimension(n) :: y
148 real(kind=rp), intent(inout), dimension(n) :: z
149 real(kind=rp), intent(inout), dimension(n) :: u
150 real(kind=rp), intent(inout), dimension(n) :: v
151 real(kind=rp), intent(inout), dimension(n) :: w
152 type(time_state_t), intent(in), optional :: time
153 integer :: i, m, k, idx(4), facet
154 real(kind=rp) :: normal(3), area
155
156 m = this%unique_mask(0)
157
158 do i = 1, m
159 k = this%unique_mask(i)
160 x(k) = u(k) * this%nx%x(i)
161 y(k) = v(k) * this%ny%x(i)
162 z(k) = w(k) * this%nz%x(i)
163 end do
164
165 end subroutine facet_normal_apply_surfvec
166
168 subroutine facet_normal_apply_surfvec_dev(this, x_d, y_d, z_d, &
169 u_d, v_d, w_d, time, strm)
170 class(facet_normal_t), intent(in), target :: this
171 type(c_ptr) :: x_d, y_d, z_d, u_d, v_d, w_d
172 type(time_state_t), intent(in), optional :: time
173 type(c_ptr), optional :: strm
174 type(c_ptr) :: strm_
175 integer :: n, m
176
177 n = this%coef%dof%size()
178 m = this%unique_mask(0)
179
180 if (present(strm)) then
181 strm_ = strm
182 else
183 strm_ = glb_cmd_queue
184 end if
185
186 if (m .gt. 0) then
187 call device_masked_gather_copy_0(this%work%x_d, u_d, &
188 this%unique_mask_d, n, m, strm_)
189 call device_col2(this%work%x_d, this%nx%x_d, m, strm_)
190 call device_masked_scatter_copy_0(x_d, this%work%x_d, &
191 this%unique_mask_d, n, m, strm_)
192 call device_masked_gather_copy_0(this%work%x_d, v_d, &
193 this%unique_mask_d, n, m, strm_)
194 call device_col2(this%work%x_d, this%ny%x_d, m, strm_)
195 call device_masked_scatter_copy_0(y_d, this%work%x_d, &
196 this%unique_mask_d, n, m, strm_)
197 call device_masked_gather_copy_0(this%work%x_d, w_d, &
198 this%unique_mask_d, n, m, strm_)
199 call device_col2(this%work%x_d, this%nz%x_d, m, strm)
200 call device_masked_scatter_copy_0(z_d, this%work%x_d, &
201 this%unique_mask_d, n, m, strm_)
202 end if
203
204 end subroutine facet_normal_apply_surfvec_dev
205
207 subroutine facet_normal_free(this)
208 class(facet_normal_t), target, intent(inout) :: this
209
210 call this%free_base()
211 if (allocated(this%unique_mask)) then
212 deallocate(this%unique_mask)
213 end if
214 if (c_associated(this%unique_mask_d)) then
215 call device_free(this%unique_mask_d)
216 end if
217
218 call this%nx%free()
219 call this%ny%free()
220 call this%nz%free()
221 call this%work%free()
222
223 end subroutine facet_normal_free
224
226 subroutine facet_normal_finalize(this, only_facets)
227 class(facet_normal_t), target, intent(inout) :: this
228 logical, optional, intent(in) :: only_facets
229 logical :: only_facets_
230 type(htable_i4_t) :: unique_point_idx
231 integer :: htable_data, rcode, i, j, idx(4), facet
232 real(kind=rp) :: area, normal(3)
233
234 if (present(only_facets)) then
235 if (.not. only_facets) then
236 call neko_error("For facet_normal_t, only_facets has to be true.")
237 end if
238 end if
239
240 call this%finalize_base(.true.)
241 ! This part is purely needed to ensure that contributions
242 ! for all faces a point is on is properly summed up.
243 ! If one simply uses the original mask, if a point is on a corner
244 ! where both faces are on the boundary
245 ! one will only get the contribution from one face, not both
246 ! We solve this by adding up the normals of both faces for these points
247 ! and storing this sum in this%nx, this%ny, this%nz.
248 ! As both contrbutions are added already,
249 ! we also ensure that we only visit each point once
250 ! and create a new mask with only unique points (this%unique_mask).
251 if (allocated(this%unique_mask)) then
252 deallocate(this%unique_mask)
253 end if
254 if (c_associated(this%unique_mask_d)) then
255 call device_free(this%unique_mask_d)
256 end if
257
258 call unique_point_idx%init(this%msk(0), htable_data)
259 j = 0
260 do i = 1, this%msk(0)
261 if (unique_point_idx%get(this%msk(i), htable_data) .ne. 0) then
262 j = j + 1
263 htable_data = j
264 call unique_point_idx%set(this%msk(i), j)
265 end if
266 end do
267
268 ! Only allocate work vectors if size is non-zero
269 if (unique_point_idx%num_entries() .gt. 0 ) then
270 call this%nx%init(unique_point_idx%num_entries())
271 call this%ny%init(unique_point_idx%num_entries())
272 call this%nz%init(unique_point_idx%num_entries())
273 call this%work%init(unique_point_idx%num_entries())
274 end if
275 allocate(this%unique_mask(0:unique_point_idx%num_entries()))
276
277 this%unique_mask(0) = unique_point_idx%num_entries()
278 do i = 1, this%unique_mask(0)
279 this%unique_mask(i) = 0
280 end do
281
282
283 do i = 1, this%msk(0)
284 rcode = unique_point_idx%get(this%msk(i), htable_data)
285 if (rcode .ne. 0) call neko_error("Facet normal: htable get failed.")
286 this%unique_mask(htable_data) = this%msk(i)
287 facet = this%facet(i)
288
289 idx = nonlinear_index(this%msk(i), this%Xh%lx, this%Xh%lx, this%Xh%lx)
290 normal = this%coef%get_normal(idx(1), idx(2), idx(3), idx(4), facet)
291 area = this%coef%get_area(idx(1), idx(2), idx(3), idx(4), facet)
292 normal = normal * area !Scale normal by area
293 this%nx%x(htable_data) = this%nx%x(htable_data) + normal(1)
294 this%ny%x(htable_data) = this%ny%x(htable_data) + normal(2)
295 this%nz%x(htable_data) = this%nz%x(htable_data) + normal(3)
296 end do
297
298 if (neko_bcknd_device .eq. 1 .and. &
299 (unique_point_idx%num_entries() .gt. 0 )) then
300 call device_map(this%unique_mask, this%unique_mask_d, &
301 size(this%unique_mask))
302 call device_memcpy(this%unique_mask, this%unique_mask_d, &
303 size(this%unique_mask), host_to_device, sync = .true.)
304 call device_memcpy(this%nx%x, this%nx%x_d, &
305 this%nx%size(), host_to_device, sync = .true.)
306 call device_memcpy(this%ny%x, this%ny%x_d, &
307 this%ny%size(), host_to_device, sync = .true.)
308 call device_memcpy(this%nz%x, this%nz%x_d, &
309 this%nz%size(), host_to_device, sync = .true.)
310 end if
311
312 call unique_point_idx%free()
313
314 end subroutine facet_normal_finalize
315
316end module facet_normal
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
subroutine, public device_masked_scatter_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Scatter a masked vector .
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:219
integer, parameter, public device_to_host
Definition device.F90:47
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:51
Dirichlet condition applied in the facet normal direction.
subroutine facet_normal_init_from_components(this, coef)
Constructor from components.
subroutine facet_normal_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
No-op vector apply on device.
subroutine facet_normal_finalize(this, only_facets)
Finalize.
subroutine facet_normal_apply_scalar_dev(this, x_d, time, strong, strm)
No-op scalar apply on device.
subroutine facet_normal_init(this, coef, json)
Constructor.
subroutine facet_normal_apply_scalar(this, x, n, time, strong)
No-op scalar apply.
subroutine facet_normal_apply_surfvec_dev(this, x_d, y_d, z_d, u_d, v_d, w_d, time, strm)
Apply in facet normal direction (vector valued, device version)
subroutine facet_normal_apply_surfvec(this, x, y, z, u, v, w, n, time)
Apply in facet normal direction (vector valued)
subroutine facet_normal_apply_vector(this, x, y, z, n, time, strong)
No-op vector apply.
subroutine facet_normal_free(this)
Destructor.
Implements a hash table ADT.
Definition htable.f90:36
Definition math.f90:60
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:397
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
Defines a vector.
Definition vector.f90:34
Base type for a boundary condition.
Definition bc.f90:62
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
Dirichlet condition in facet normal direction.
Integer based hash table.
Definition htable.f90:82
A struct that contains all info about the time, expand as needed.