44 use json_module,
only : json_file
45 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
55 integer,
allocatable :: unique_mask(:)
56 type(c_ptr) :: unique_mask_d = c_null_ptr
64 procedure, pass(this) :: apply_surfvec_dev => &
69 procedure, pass(this) :: init_from_components => &
84 type(
coef_t),
target,
intent(in) :: coef
85 type(json_file),
intent(inout) :: json
87 call this%init_from_components(coef)
94 type(
coef_t),
target,
intent(in) :: coef
96 call this%init_base(coef)
102 integer,
intent(in) :: n
103 real(kind=
rp),
intent(inout),
dimension(n) :: x
105 logical,
intent(in),
optional :: strong
111 type(c_ptr),
intent(inout) :: x_d
113 logical,
intent(in),
optional :: strong
114 type(c_ptr),
intent(inout) :: strm
122 type(c_ptr),
intent(inout) :: x_d
123 type(c_ptr),
intent(inout) :: y_d
124 type(c_ptr),
intent(inout) :: z_d
126 logical,
intent(in),
optional :: strong
127 type(c_ptr),
intent(inout) :: strm
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
139 logical,
intent(in),
optional :: strong
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
153 integer :: i, m, k, idx(4), facet
154 real(kind=
rp) :: normal(3), area
156 m = this%unique_mask(0)
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)
169 u_d, v_d, w_d, time, strm)
171 type(c_ptr) :: x_d, y_d, z_d, u_d, v_d, w_d
173 type(c_ptr),
optional :: strm
177 n = this%coef%dof%size()
178 m = this%unique_mask(0)
180 if (
present(strm))
then
188 this%unique_mask_d, n, m, strm_)
189 call device_col2(this%work%x_d, this%nx%x_d, m, strm_)
191 this%unique_mask_d, n, m, strm_)
193 this%unique_mask_d, n, m, strm_)
194 call device_col2(this%work%x_d, this%ny%x_d, m, strm_)
196 this%unique_mask_d, n, m, strm_)
198 this%unique_mask_d, n, m, strm_)
199 call device_col2(this%work%x_d, this%nz%x_d, m, strm)
201 this%unique_mask_d, n, m, strm_)
210 call this%free_base()
211 if (
allocated(this%unique_mask))
then
212 deallocate(this%unique_mask)
214 if (c_associated(this%unique_mask_d))
then
221 call this%work%free()
228 logical,
optional,
intent(in) :: only_facets
229 logical :: only_facets_
231 integer :: htable_data, rcode, i, j, idx(4), facet
232 real(kind=
rp) :: area, normal(3)
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.")
240 call this%finalize_base(.true.)
251 if (
allocated(this%unique_mask))
then
252 deallocate(this%unique_mask)
254 if (c_associated(this%unique_mask_d))
then
258 call unique_point_idx%init(this%msk(0), htable_data)
260 do i = 1, this%msk(0)
261 if (unique_point_idx%get(this%msk(i), htable_data) .ne. 0)
then
264 call unique_point_idx%set(this%msk(i), j)
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())
275 allocate(this%unique_mask(0:unique_point_idx%num_entries()))
277 this%unique_mask(0) = unique_point_idx%num_entries()
278 do i = 1, this%unique_mask(0)
279 this%unique_mask(i) = 0
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)
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
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)
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))
312 call unique_point_idx%free()
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Map a Fortran array to a device (allocate and associate)
Copy data between host and device (or device and device)
Defines a boundary condition.
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.
integer, parameter, public host_to_device
subroutine, public device_free(x_d)
Deallocate memory on the device.
integer, parameter, public device_to_host
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
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.
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Module with things related to the simulation time.
Base type for a boundary condition.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Dirichlet condition in facet normal direction.
Integer based hash table.
A struct that contains all info about the time, expand as needed.