39 use iso_c_binding,
only : c_associated
52 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
53 use json_module,
only : json_file
62 type,
public,
abstract ::
bc_t
64 integer,
allocatable :: msk(:)
66 integer,
allocatable :: facet(:)
70 type(
coef_t),
pointer :: coef => null()
72 type(
mesh_t),
pointer :: msh => null()
78 type(c_ptr) :: msk_d = c_null_ptr
80 type(c_ptr) :: facet_d = c_null_ptr
85 logical :: strong = .true.
88 logical :: updated = .false.
90 character(len=:),
allocatable :: name
92 integer,
allocatable :: zone_indices(:)
134 class(
bc_t),
pointer :: ptr => null()
139 class(
bc_t),
allocatable :: obj
147 class(
bc_t),
intent(inout),
target :: this
148 type(
coef_t),
target,
intent(in) :: coef
149 type(json_file),
intent(inout) :: json
157 class(
bc_t),
intent(inout),
target :: this
165 class(
bc_t),
intent(inout),
target :: this
166 logical,
optional,
intent(in) :: only_facets
179 class(
bc_t),
intent(inout) :: this
180 integer,
intent(in) :: n
181 real(kind=
rp),
intent(inout),
dimension(n) :: x
183 logical,
intent(in),
optional :: strong
199 class(
bc_t),
intent(inout) :: this
200 integer,
intent(in) :: n
201 real(kind=
rp),
intent(inout),
dimension(n) :: x
202 real(kind=
rp),
intent(inout),
dimension(n) :: y
203 real(kind=
rp),
intent(inout),
dimension(n) :: z
205 logical,
intent(in),
optional :: strong
219 class(
bc_t),
intent(inout),
target :: this
220 type(c_ptr),
intent(inout) :: x_d
222 logical,
intent(in),
optional :: strong
223 type(c_ptr),
intent(inout) :: strm
238 class(
bc_t),
intent(inout),
target :: this
239 type(c_ptr),
intent(inout) :: x_d
240 type(c_ptr),
intent(inout) :: y_d
241 type(c_ptr),
intent(inout) :: z_d
243 logical,
intent(in),
optional :: strong
244 type(c_ptr),
intent(inout) :: strm
253 class(
bc_t),
intent(inout) :: this
254 type(
coef_t),
target,
intent(in) :: coef
260 this%Xh => this%dof%Xh
261 this%msh => this%dof%msh
263 call this%marked_facet%init()
269 class(
bc_t),
intent(inout) :: this
271 call this%marked_facet%free()
278 if (
allocated(this%msk))
then
285 if (
allocated(this%facet))
then
289 deallocate(this%facet)
292 if (
allocated(this%name))
then
293 deallocate(this%name)
296 if (
allocated(this%zone_indices))
then
297 deallocate(this%zone_indices)
311 class(
bc_t),
intent(inout) :: this
312 type(
field_t),
intent(inout) :: x
313 type(
field_t),
intent(inout) :: y
314 type(
field_t),
intent(inout) :: z
316 logical,
intent(in),
optional :: strong
317 type(c_ptr),
intent(inout),
optional :: strm
320 character(len=256) :: msg
326 if (y%size() .ne. n .or. z%size() .ne. n)
then
327 msg =
"Fields x, y, z must have the same size in " // &
328 "bc_list_apply_vector_field"
334 if (
present(strm))
then
340 call this%apply_vector_dev(x%x_d, y%x_d, z%x_d, time = time, &
341 strong = strong, strm = strm_)
343 call this%apply_vector(x%x, y%x, z%x, n, time = time, strong = strong)
355 class(
bc_t),
intent(inout) :: this
356 type(
field_t),
intent(inout) :: x
358 logical,
intent(in),
optional :: strong
359 type(c_ptr),
intent(inout),
optional :: strm
368 if (
present(strm))
then
374 call this%apply_scalar_dev(x%x_d, time = time, strong = strong, &
377 call this%apply_scalar(x%x, n, time = time)
386 class(
bc_t),
intent(inout) :: this
387 integer,
intent(in) :: facet
388 integer,
intent(in) :: el
392 call this%marked_facet%push(t)
399 class(
bc_t),
intent(inout) :: this
404 fp => facet_list%array()
405 do i = 1, facet_list%size()
406 call this%marked_facet%push(fp(i))
414 class(
bc_t),
intent(inout) :: this
417 do i = 1, bc_zone%size
418 call this%marked_facet%push(bc_zone%facet_el(i))
430 class(
bc_t),
target,
intent(inout) :: this
431 logical,
optional,
intent(in) :: only_facets
435 integer :: facet_size, facet, el
436 logical :: only_facet = .false.
437 integer :: i, j, k, l, msk_c
438 integer :: lx, ly, lz, n
439 character(len=LOG_SIZE) :: log_buf
443 if (
present(only_facets))
then
444 only_facet = only_facets
452 allocate(this%msk(0:facet_size * this%marked_facet%size()))
453 allocate(this%facet(0:facet_size * this%marked_facet%size()))
456 bfp => this%marked_facet%array()
462 do i = 1, this%marked_facet%size()
464 facet = bc_facet%x(1)
472 this%facet(msk_c) = 1
479 this%msk(msk_c) =
linear_index(lx, k, l, el, lx, ly, lz)
480 this%facet(msk_c) = 2
488 this%facet(msk_c) = 3
495 this%msk(msk_c) =
linear_index(j, ly, l, el, lx, ly, lz)
496 this%facet(msk_c) = 4
504 this%facet(msk_c) = 5
511 this%msk(msk_c) =
linear_index(j, k, lz, el, lx, ly, lz)
512 this%facet(msk_c) = 6
517 this%facet(0) = msk_c
525 if ( .not. only_facet)
then
527 call test_field%init(this%dof)
529 n = test_field%size()
530 test_field%x = 0.0_rp
533 test_field%x(this%msk(i),1,1,1) = 1.0
540 call this%coef%gs_h%op(test_field,
gs_op_add)
546 do i = 1, this%dof%size()
547 if (test_field%x(i,1,1,1) .gt. 0.5)
then
553 allocate(this%msk(0:msk_c))
555 do i = 1, this%dof%size()
556 if (test_field%x(i,1,1,1) .gt. 0.5)
then
562 call test_field%free()
574 if (.not.
allocated(this%name))
then
578 write(log_buf,
'(A,A)')
'BC assigned name : ', trim(this%name)
594 class(
bc_t),
intent(inout) :: this
595 character(len=*),
intent(in) :: file_name
600 call bdry_field%init(this%coef%dof,
'bdry')
604 bdry_field%x(k,1,1,1) = 1.0_rp
606 call dump_file%init(file_name)
607 call dump_file%write(bdry_field)
Apply the boundary condition to a scalar field on the device.
Apply the boundary condition to a scalar field.
Apply the boundary condition to a vector field on the device.
Apply the boundary condition to a vector field.
Finalize by building the mask and facet arrays.
Map a Fortran array to a device (allocate and associate)
Copy data between host and device (or device and device)
Unmap a Fortran array from a device (deassociate and free)
Defines a boundary condition.
subroutine bc_mark_zone(this, bc_zone)
Mark all facets from a zone.
subroutine bc_free_base(this)
Destructor for the base type, bc_t.
subroutine bc_init_base(this, coef)
Constructor.
subroutine bc_apply_scalar_generic(this, x, time, strong, strm)
Apply the boundary condition to a scalar field. Dispatches to the CPU or the device version.
subroutine bc_finalize_base(this, only_facets)
Finalize the construction of the bc by populting the msk and facet arrays.
subroutine bc_apply_vector_generic(this, x, y, z, time, strong, strm)
Apply the boundary condition to a vector field. Dispatches to the CPU or the device version.
subroutine bc_debug_mask(this, file_name)
Write a field showing the mask of the bc.
subroutine bc_mark_facet(this, facet, el)
Mark facet on element el as part of the boundary condition.
subroutine bc_mark_facets(this, facet_list)
Mark all facets from a (facet, el) tuple list.
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
integer, parameter, public device_to_host
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Defines a mapping of the degrees of freedom.
Defines a zone as a subset of facets in a mesh.
Module for file I/O operations.
Defines Gather-scatter operations.
integer, parameter, public gs_op_add
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Defines a function space.
Implements a dynamic stack ADT.
Module with things related to the simulation time.
character(len=100) function, dimension(:), allocatable, public split_string(string, delimiter)
Split a string based on delimiter (tokenizer) OBS: very hacky, this should really be improved,...
pure integer function, public linear_index(i, j, k, l, lx, ly, lz)
Compute the address of a (i,j,k,l) array with sizes (1:lx, 1:ly, 1:lz, :)
Base type for a boundary condition.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
The function space for the SEM solution fields.
Integer 2-tuple based stack.
A struct that contains all info about the time, expand as needed.