46 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
51 type,
public,
abstract ::
bc_t
53 integer,
allocatable :: msk(:)
55 integer,
allocatable :: facet(:)
67 type(c_ptr) :: msk_d = c_null_ptr
69 type(c_ptr) :: facet_d = c_null_ptr
105 class(
bc_t),
allocatable :: obj
117 class(
bc_t),
intent(inout) :: this
118 integer,
intent(in) :: n
119 real(kind=
rp),
intent(inout),
dimension(n) :: x
120 real(kind=
rp),
intent(in),
optional :: t
121 integer,
intent(in),
optional :: tstep
136 class(
bc_t),
intent(inout) :: this
137 integer,
intent(in) :: n
138 real(kind=
rp),
intent(inout),
dimension(n) :: x
139 real(kind=
rp),
intent(inout),
dimension(n) :: y
140 real(kind=
rp),
intent(inout),
dimension(n) :: z
141 real(kind=
rp),
intent(in),
optional :: t
142 integer,
intent(in),
optional :: tstep
150 class(
bc_t),
intent(inout),
target :: this
161 class(
bc_t),
intent(inout),
target :: this
163 real(kind=
rp),
intent(in),
optional :: t
164 integer,
intent(in),
optional :: tstep
177 class(
bc_t),
intent(inout),
target :: this
181 real(kind=
rp),
intent(in),
optional :: t
182 integer,
intent(in),
optional :: tstep
191 class(
bc_t),
intent(inout) :: this
192 type(
coef_t),
target,
intent(in) :: coef
198 this%Xh => this%dof%Xh
199 this%msh => this%dof%msh
201 call this%marked_facet%init()
207 class(
bc_t),
intent(inout) :: this
209 call this%marked_facet%free()
216 if (
allocated(this%msk))
then
220 if (
allocated(this%facet))
then
221 deallocate(this%facet)
224 if (c_associated(this%msk_d))
then
226 this%msk_d = c_null_ptr
229 if (c_associated(this%facet_d))
then
231 this%facet_d = c_null_ptr
240 class(
bc_t),
intent(inout) :: this
241 integer,
intent(in) :: facet
242 integer,
intent(in) :: el
246 call this%marked_facet%push(t)
253 class(
bc_t),
intent(inout) :: this
258 fp => facet_list%array()
259 do i = 1, facet_list%size()
260 call this%marked_facet%push(fp(i))
268 class(
bc_t),
intent(inout) :: this
271 do i = 1, bc_zone%size
272 call this%marked_facet%push(bc_zone%facet_el(i))
283 class(
bc_t),
intent(inout) :: this
285 character(len=*) :: bc_key
286 character(len=100),
allocatable :: split_key(:)
287 character(len=NEKO_MSH_MAX_ZLBL_LEN) :: bc_labels(:)
288 integer :: i, j, k, l, msh_bc_type
291 if(trim(bc_key) .eq.
'o' .or. trim(bc_key) .eq.
'on' &
292 .or. trim(bc_key) .eq.
'o+dong' .or. trim(bc_key) .eq.
'on+dong')
then
294 else if(trim(bc_key) .eq.
'd_pres')
then
296 else if(trim(bc_key) .eq.
'w')
then
298 else if(trim(bc_key) .eq.
'v')
then
300 else if(trim(bc_key) .eq.
'd_vel_u')
then
302 else if(trim(bc_key) .eq.
'd_vel_v')
then
304 else if(trim(bc_key) .eq.
'd_vel_w')
then
306 else if(trim(bc_key) .eq.
'sym')
then
308 else if(trim(bc_key) .eq.
'sh')
then
310 else if(trim(bc_key) .eq.
'wm')
then
314 do i = 1,
size(bc_labels)
317 if (index(trim(bc_labels(i)),
'/') .eq. 0)
then
318 if (trim(bc_key) .eq. trim(bc_labels(i)))
then
321 do j = 1,this%msh%nelv
322 do k = 1, 2 * this%msh%gdim
323 if (this%msh%facet_type(k,j) .eq. -i)
then
324 this%msh%facet_type(k,j) = msh_bc_type
331 do l = 1,
size(split_key)
332 if (trim(split_key(l)) .eq. trim(bc_key))
then
335 do j = 1,this%msh%nelv
336 do k = 1, 2 * this%msh%gdim
337 if (this%msh%facet_type(k,j) .eq. -i)
then
338 this%msh%facet_type(k,j) = msh_bc_type
353 class(
bc_t),
target,
intent(inout) :: this
356 integer :: facet_size, facet, el
357 integer :: i, j, k, l, msk_c
358 integer :: lx, ly, lz, n
368 allocate(this%msk(0:facet_size * this%marked_facet%size()))
369 allocate(this%facet(0:facet_size * this%marked_facet%size()))
372 bfp => this%marked_facet%array()
378 do i = 1, this%marked_facet%size()
380 facet = bc_facet%x(1)
388 this%facet(msk_c) = 1
396 this%facet(msk_c) = 2
404 this%facet(msk_c) = 3
412 this%facet(msk_c) = 4
420 this%facet(msk_c) = 5
428 this%facet(msk_c) = 6
435 this%facet(0) = msk_c
438 n = facet_size * this%marked_facet%size() + 1
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.
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_finalize(this)
Finalize the construction of the bc by populting the msk and facet arrays.
subroutine bc_free_base(this)
Destructor for the base type, bc_t.
subroutine bc_init_base(this, coef)
Constructor.
subroutine bc_mark_zones_from_list(this, bc_zones, bc_key, bc_labels)
Mark all facets from a list of zones, also marks type of bc in the mesh. The facet_type in mesh is be...
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
subroutine, public device_free(x_d)
Deallocate memory on the device.
Defines a mapping of the degrees of freedom.
Defines a zone as a subset of facets in a mesh.
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, :)
Wrapper around an allocatable bc.
Pointer to boundary condtiion.
Base type for a boundary condition.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
The function space for the SEM solution fields.
Integer 2-tuple based stack.