36 use num_types,
only : rp
37 use neko_config,
only : neko_bcknd_device
38 use vector,
only : vector_t
39 use coefs,
only : coef_t
41 use utils,
only : neko_error, nonlinear_index
42 use json_module,
only : json_file
43 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
44 use htable,
only : htable_i4_t
45 use device,
only : device_map, device_memcpy, device_free, &
47 use time_state,
only : time_state_t
53 integer,
allocatable :: unique_mask(:)
54 type(c_ptr) :: unique_mask_d = c_null_ptr
55 type(vector_t) :: nx, ny, nz, work
57 procedure, pass(this) :: apply_scalar => normal_vec_bcs_apply_scalar
58 procedure, pass(this) :: apply_scalar_dev => &
59 normal_vec_bcs_apply_scalar_dev
60 procedure, pass(this) :: apply_vector => normal_vec_bcs_apply_vector
61 procedure, pass(this) :: apply_vector_dev => &
62 normal_vec_bcs_apply_vector_dev
63 procedure, pass(this) :: apply_n_dot => normal_vec_bcs_apply_n_dot
64 procedure, pass(this) :: apply_n_cross => normal_vec_bcs_apply_n_cross
66 procedure, pass(this) :: init => normal_vec_bcs_init
68 procedure, pass(this) :: init_from_components => &
69 normal_vec_bcs_init_from_components
71 procedure, pass(this) :: free => normal_vec_bcs_free
73 procedure, pass(this) :: finalize => normal_vec_bcs_finalize
82 subroutine normal_vec_bcs_init(this, coef, json)
84 type(coef_t),
target,
intent(in) :: coef
85 type(json_file),
intent(inout) ::json
87 call this%init_from_components(coef)
88 end subroutine normal_vec_bcs_init
93 subroutine normal_vec_bcs_init_from_components(this, coef)
95 type(coef_t),
target,
intent(in) :: coef
97 call this%init_base(coef)
98 end subroutine normal_vec_bcs_init_from_components
106 subroutine normal_vec_bcs_apply_scalar(this, x, n, time, strong)
108 integer,
intent(in) :: n
109 real(kind=rp),
intent(inout),
dimension(n) :: x
110 type(time_state_t),
intent(in),
optional :: time
111 logical,
intent(in),
optional :: strong
112 end subroutine normal_vec_bcs_apply_scalar
120 subroutine normal_vec_bcs_apply_scalar_dev(this, x_d, time, strong, strm)
122 type(c_ptr),
intent(inout) :: x_d
123 type(time_state_t),
intent(in),
optional :: time
124 logical,
intent(in),
optional :: strong
125 type(c_ptr),
intent(inout) :: strm
127 end subroutine normal_vec_bcs_apply_scalar_dev
137 subroutine normal_vec_bcs_apply_vector_dev(this, x_d, y_d, z_d, time, &
140 type(c_ptr),
intent(inout) :: x_d
141 type(c_ptr),
intent(inout) :: y_d
142 type(c_ptr),
intent(inout) :: z_d
143 type(time_state_t),
intent(in),
optional :: time
144 logical,
intent(in),
optional :: strong
145 type(c_ptr),
intent(inout) :: strm
147 end subroutine normal_vec_bcs_apply_vector_dev
157 subroutine normal_vec_bcs_apply_vector(this, x, y, z, n, time, strong)
159 integer,
intent(in) :: n
160 real(kind=rp),
intent(inout),
dimension(n) :: x
161 real(kind=rp),
intent(inout),
dimension(n) :: y
162 real(kind=rp),
intent(inout),
dimension(n) :: z
163 type(time_state_t),
intent(in),
optional :: time
164 logical,
intent(in),
optional :: strong
165 end subroutine normal_vec_bcs_apply_vector
175 subroutine normal_vec_bcs_apply_n_dot(this, f, u, v, w, n, time)
177 integer,
intent(in) :: n
178 real(kind=rp),
intent(inout),
dimension(n) :: f
179 real(kind=rp),
intent(inout),
dimension(n) :: u
180 real(kind=rp),
intent(inout),
dimension(n) :: v
181 real(kind=rp),
intent(inout),
dimension(n) :: w
182 type(time_state_t),
intent(in),
optional :: time
185 m = this%unique_mask(0)
190 k = this%unique_mask(i)
191 f(k) = u(k) * this%nx%x(i) + v(k) * this%ny%x(i) + w(k) * this%nz%x(i)
194 end subroutine normal_vec_bcs_apply_n_dot
206 subroutine normal_vec_bcs_apply_n_cross(this, x, y, z, u, v, w, n, time)
208 integer,
intent(in) :: n
209 real(kind=rp),
intent(inout),
dimension(n) :: x
210 real(kind=rp),
intent(inout),
dimension(n) :: y
211 real(kind=rp),
intent(inout),
dimension(n) :: z
212 real(kind=rp),
intent(inout),
dimension(n) :: u
213 real(kind=rp),
intent(inout),
dimension(n) :: v
214 real(kind=rp),
intent(inout),
dimension(n) :: w
215 type(time_state_t),
intent(in),
optional :: time
218 m = this%unique_mask(0)
221 k = this%unique_mask(i)
222 x(k) = this%ny%x(i) * w(k) - this%nz%x(i) * v(k)
223 y(k) = this%nz%x(i) * u(k) - this%nx%x(i) * w(k)
224 z(k) = this%nx%x(i) * v(k) - this%ny%x(i) * u(k)
227 end subroutine normal_vec_bcs_apply_n_cross
231 subroutine normal_vec_bcs_free(this)
234 call this%free_base()
235 if (
allocated(this%unique_mask))
then
236 deallocate(this%unique_mask)
238 if (c_associated(this%unique_mask_d))
then
239 call device_free(this%unique_mask_d)
245 call this%work%free()
247 end subroutine normal_vec_bcs_free
252 subroutine normal_vec_bcs_finalize(this, only_facets)
254 logical,
optional,
intent(in) :: only_facets
255 type(htable_i4_t) :: unique_point_idx
256 integer :: htable_data, rcode, i, j, idx(4), facet
257 real(kind=rp) :: area, normal(3)
259 if (
present(only_facets))
then
260 if (only_facets .eqv. .false.)
then
261 call neko_error(
"For normal_vec_bcs_t, only_facets has to be true.")
265 call this%finalize_base(.true.)
276 if (
allocated(this%unique_mask))
then
277 deallocate(this%unique_mask)
279 if (c_associated(this%unique_mask_d))
then
280 call device_free(this%unique_mask_d)
283 call unique_point_idx%init(this%msk(0), htable_data)
285 do i = 1, this%msk(0)
286 if (unique_point_idx%get(this%msk(i),htable_data) .ne. 0)
then
289 call unique_point_idx%set(this%msk(i), j)
294 if (unique_point_idx%num_entries() .gt. 0 )
then
295 call this%nx%init(unique_point_idx%num_entries())
296 call this%ny%init(unique_point_idx%num_entries())
297 call this%nz%init(unique_point_idx%num_entries())
298 call this%work%init(unique_point_idx%num_entries())
300 allocate(this%unique_mask(0:unique_point_idx%num_entries()))
302 this%unique_mask(0) = unique_point_idx%num_entries()
303 do i = 1, this%unique_mask(0)
304 this%unique_mask(i) = 0
308 do i = 1, this%msk(0)
309 rcode = unique_point_idx%get(this%msk(i), htable_data)
310 if (rcode .ne. 0)
call neko_error(
"Facet normal: htable get failed.")
311 this%unique_mask(htable_data) = this%msk(i)
312 facet = this%facet(i)
314 idx = nonlinear_index(this%msk(i), this%Xh%lx, this%Xh%lx, this%Xh%lx)
315 normal = this%coef%get_normal(idx(1), idx(2), idx(3), idx(4), facet)
316 area = this%coef%get_area(idx(1), idx(2), idx(3), idx(4), facet)
317 normal = normal * area
318 this%nx%x(htable_data) = this%nx%x(htable_data) + normal(1)
319 this%ny%x(htable_data) = this%ny%x(htable_data) + normal(2)
320 this%nz%x(htable_data) = this%nz%x(htable_data) + normal(3)
323 if (neko_bcknd_device .eq. 1 .and. &
324 (unique_point_idx%num_entries() .gt. 0 ))
then
325 call device_map(this%unique_mask, this%unique_mask_d, &
326 size(this%unique_mask))
327 call device_memcpy(this%unique_mask, this%unique_mask_d, &
328 size(this%unique_mask), host_to_device, sync = .true.)
329 call device_memcpy(this%nx%x, this%nx%x_d, &
330 this%nx%size(), host_to_device, sync = .true.)
331 call device_memcpy(this%ny%x, this%ny%x_d, &
332 this%ny%size(), host_to_device, sync = .true.)
333 call device_memcpy(this%nz%x, this%nz%x_d, &
334 this%nz%size(), host_to_device, sync = .true.)
337 call unique_point_idx%free()
339 end subroutine normal_vec_bcs_finalize
341end module normal_vec_bcs
Dirichlet condition in facet normal direction.