39 use iso_c_binding,
only: c_associated
51 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
52 use json_module,
only : json_file
57 type,
public,
abstract ::
bc_t
59 integer,
allocatable :: msk(:)
61 integer,
allocatable :: facet(:)
65 type(
coef_t),
pointer :: coef => null()
67 type(
mesh_t),
pointer :: msh => null()
73 type(c_ptr) :: msk_d = c_null_ptr
75 type(c_ptr) :: facet_d = c_null_ptr
80 logical :: strong = .true.
57 type,
public,
abstract ::
bc_t …
122 class(
bc_t),
pointer :: ptr => null()
127 class(
bc_t),
allocatable :: obj
135 class(
bc_t),
intent(inout),
target :: this
136 type(
coef_t),
intent(in) :: coef
137 type(json_file),
intent(inout) :: json
145 class(
bc_t),
intent(inout),
target :: this
153 class(
bc_t),
intent(inout),
target :: this
167 class(
bc_t),
intent(inout) :: this
168 integer,
intent(in) :: n
169 real(kind=
rp),
intent(inout),
dimension(n) :: x
170 real(kind=
rp),
intent(in),
optional :: t
171 integer,
intent(in),
optional :: tstep
172 logical,
intent(in),
optional :: strong
188 class(
bc_t),
intent(inout) :: this
189 integer,
intent(in) :: n
190 real(kind=
rp),
intent(inout),
dimension(n) :: x
191 real(kind=
rp),
intent(inout),
dimension(n) :: y
192 real(kind=
rp),
intent(inout),
dimension(n) :: z
193 real(kind=
rp),
intent(in),
optional :: t
194 integer,
intent(in),
optional :: tstep
195 logical,
intent(in),
optional :: strong
209 class(
bc_t),
intent(inout),
target :: this
211 real(kind=
rp),
intent(in),
optional :: t
212 integer,
intent(in),
optional :: tstep
213 logical,
intent(in),
optional :: strong
229 class(
bc_t),
intent(inout),
target :: this
233 real(kind=
rp),
intent(in),
optional :: t
234 integer,
intent(in),
optional :: tstep
235 logical,
intent(in),
optional :: strong
244 class(
bc_t),
intent(inout) :: this
245 type(
coef_t),
target,
intent(in) :: coef
251 this%Xh => this%dof%Xh
252 this%msh => this%dof%msh
254 call this%marked_facet%init()
260 class(
bc_t),
intent(inout) :: this
262 call this%marked_facet%free()
269 if (
allocated(this%msk))
then
273 if (
allocated(this%facet))
then
274 deallocate(this%facet)
277 if (c_associated(this%msk_d))
then
279 this%msk_d = c_null_ptr
282 if (c_associated(this%facet_d))
then
284 this%facet_d = c_null_ptr
298 class(
bc_t),
intent(inout) :: this
299 integer,
intent(in) :: n
300 real(kind=
rp),
intent(inout),
dimension(n) :: x
301 real(kind=
rp),
intent(inout),
dimension(n) :: y
302 real(kind=
rp),
intent(inout),
dimension(n) :: z
303 real(kind=
rp),
intent(in),
optional :: t
304 integer,
intent(in),
optional :: tstep
315 if (
present(t) .and.
present(tstep))
then
316 call this%apply_vector_dev(x_d, y_d, z_d, t = t, tstep = tstep)
317 else if (
present(t))
then
318 call this%apply_vector_dev(x_d, y_d, z_d, t = t)
319 else if (
present(tstep))
then
320 call this%apply_vector_dev(x_d, y_d, z_d, tstep = tstep)
322 call this%apply_vector_dev(x_d, y_d, z_d)
325 if (
present(t) .and.
present(tstep))
then
326 call this%apply_vector(x, y, z, n, t = t, tstep = tstep)
327 else if (
present(t))
then
328 call this%apply_vector(x, y, z, n, t = t)
329 else if (
present(tstep))
then
330 call this%apply_vector(x, y, z, n, tstep = tstep)
332 call this%apply_vector(x, y, z, n)
347 class(
bc_t),
intent(inout) :: this
348 integer,
intent(in) :: n
349 real(kind=
rp),
intent(inout),
dimension(n) :: x
350 real(kind=
rp),
intent(in),
optional :: t
351 integer,
intent(in),
optional :: tstep
358 if (
present(t) .and.
present(tstep))
then
359 call this%apply_scalar_dev(x_d, t = t, tstep = tstep)
360 else if (
present(t))
then
361 call this%apply_scalar_dev(x_d, t = t)
362 else if (
present(tstep))
then
363 call this%apply_scalar_dev(x_d, tstep = tstep)
365 call this%apply_scalar_dev(x_d)
368 if (
present(t) .and.
present(tstep))
then
369 call this%apply_scalar(x, n, t = t, tstep = tstep)
370 else if (
present(t))
then
371 call this%apply_scalar(x, n, t = t)
372 else if (
present(tstep))
then
373 call this%apply_scalar(x, n, tstep = tstep)
375 call this%apply_scalar(x, n)
385 class(
bc_t),
intent(inout) :: this
386 integer,
intent(in) :: facet
387 integer,
intent(in) :: el
391 call this%marked_facet%push(t)
398 class(
bc_t),
intent(inout) :: this
403 fp => facet_list%array()
404 do i = 1, facet_list%size()
405 call this%marked_facet%push(fp(i))
413 class(
bc_t),
intent(inout) :: this
416 do i = 1, bc_zone%size
417 call this%marked_facet%push(bc_zone%facet_el(i))
429 class(
bc_t),
target,
intent(inout) :: this
430 logical,
optional,
intent(in) :: only_facets
434 integer :: facet_size, facet, el
435 logical :: only_facet = .false.
436 integer :: i, j, k, l, msk_c
437 integer :: lx, ly, lz, n
441 if (
present(only_facets))
then
442 only_facet = only_facets
450 allocate(this%msk(0:facet_size * this%marked_facet%size()))
451 allocate(this%facet(0:facet_size * this%marked_facet%size()))
454 bfp => this%marked_facet%array()
460 do i = 1, this%marked_facet%size()
462 facet = bc_facet%x(1)
470 this%facet(msk_c) = 1
477 this%msk(msk_c) =
linear_index(lx, k, l, el, lx, ly, lz)
478 this%facet(msk_c) = 2
486 this%facet(msk_c) = 3
493 this%msk(msk_c) =
linear_index(j, ly, l, el, lx, ly, lz)
494 this%facet(msk_c) = 4
502 this%facet(msk_c) = 5
509 this%msk(msk_c) =
linear_index(j, k, lz, el, lx, ly, lz)
510 this%facet(msk_c) = 6
515 if ( .not. only_facet)
then
517 call test_field%init(this%dof)
519 n = test_field%size()
520 test_field%x = 0.0_rp
523 test_field%x(this%msk(i),1,1,1) = 1.0
530 call this%coef%gs_h%op(test_field,gs_op_add)
536 do i = 1, this%dof%size()
537 if (test_field%x(i,1,1,1) .gt. 0.5)
then
543 allocate(this%msk(0:msk_c))
545 do i = 1, this%dof%size()
546 if (test_field%x(i,1,1,1) .gt. 0.5)
then
552 call test_field%free()
556 this%facet(0) = msk_c
578 class(
bc_t),
intent(inout) :: this
579 character(len=*),
intent(in) :: file_name
584 call bdry_field%init(this%coef%dof,
'bdry')
588 bdry_field%x(k,1,1,1) = 1.0_rp
590 dump_file =
file_t(file_name)
591 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.
Return the device pointer for an associated Fortran array.
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 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_apply_scalar_generic(this, x, n, t, tstep)
Apply the boundary condition to a scalar field. Dispatches to the CPU or the device version.
subroutine bc_init_base(this, coef)
Constructor.
subroutine bc_finalize_base(this, only_facets)
Finalize the construction of the bc by populting the msk and facet arrays.
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.
subroutine bc_apply_vector_generic(this, x, y, z, n, t, tstep)
Apply the boundary condition to a vector field. Dispatches to the CPU or the device version.
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
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.
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.
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.