Neko  0.8.1
A portable framework for high-order spectral element flow simulations
bc.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2024, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module bc
35  use neko_config
36  use num_types
37  use device
38  use dofmap, only : dofmap_t
39  use coefs, only : coef_t
40  use space, only : space_t
42  use facet_zone, only : facet_zone_t
43  use stack, only : stack_i4t2_t
44  use tuple, only : tuple_i4_t
46  use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
47  implicit none
48  private
49 
51  type, public, abstract :: bc_t
53  integer, allocatable :: msk(:)
55  integer, allocatable :: facet(:)
57  type(dofmap_t), pointer :: dof
59  type(coef_t), pointer :: coef
61  type(mesh_t), pointer :: msh
63  type(space_t), pointer :: xh
65  type(stack_i4t2_t) :: marked_facet
67  type(c_ptr) :: msk_d = c_null_ptr
69  type(c_ptr) :: facet_d = c_null_ptr
70  contains
72  procedure, pass(this) :: init => bc_init
74  procedure, pass(this) :: free => bc_free
76  procedure, pass(this) :: mark_facet => bc_mark_facet
78  procedure, pass(this) :: mark_facets => bc_mark_facets
80  procedure, pass(this) :: mark_zones_from_list => bc_mark_zones_from_list
82  procedure, pass(this) :: mark_zone => bc_mark_zone
85  procedure, pass(this) :: finalize => bc_finalize
87  procedure(bc_apply_scalar), pass(this), deferred :: apply_scalar
89  procedure(bc_apply_vector), pass(this), deferred :: apply_vector
91  procedure(bc_apply_scalar_dev), pass(this), deferred :: apply_scalar_dev
93  procedure(bc_apply_vector_dev), pass(this), deferred :: apply_vector_dev
94  end type bc_t
95 
97  type, private :: bcp_t
98  class(bc_t), pointer :: bcp
99  end type bcp_t
100 
102  type, public :: bc_list_t
103  type(bcp_t), allocatable :: bc(:)
105  integer :: n
107  integer :: size
108  end type bc_list_t
109 
110  abstract interface
111 
116  subroutine bc_apply_scalar(this, x, n, t, tstep)
117  import :: bc_t
118  import :: rp
119  class(bc_t), intent(inout) :: this
120  integer, intent(in) :: n
121  real(kind=rp), intent(inout), dimension(n) :: x
122  real(kind=rp), intent(in), optional :: t
123  integer, intent(in), optional :: tstep
124  end subroutine bc_apply_scalar
125  end interface
126 
127  abstract interface
128 
135  subroutine bc_apply_vector(this, x, y, z, n, t, tstep)
136  import :: bc_t
137  import :: rp
138  class(bc_t), intent(inout) :: this
139  integer, intent(in) :: n
140  real(kind=rp), intent(inout), dimension(n) :: x
141  real(kind=rp), intent(inout), dimension(n) :: y
142  real(kind=rp), intent(inout), dimension(n) :: z
143  real(kind=rp), intent(in), optional :: t
144  integer, intent(in), optional :: tstep
145  end subroutine bc_apply_vector
146  end interface
147 
148  abstract interface
149 
151  subroutine bc_apply_scalar_dev(this, x_d, t, tstep)
152  import :: c_ptr
153  import :: bc_t
154  import :: rp
155  class(bc_t), intent(inout), target :: this
156  type(c_ptr) :: x_d
157  real(kind=rp), intent(in), optional :: t
158  integer, intent(in), optional :: tstep
159  end subroutine bc_apply_scalar_dev
160  end interface
161 
162  abstract interface
163 
167  subroutine bc_apply_vector_dev(this, x_d, y_d, z_d, t, tstep)
168  import :: c_ptr
169  import :: bc_t
170  import :: rp
171  class(bc_t), intent(inout), target :: this
172  type(c_ptr) :: x_d
173  type(c_ptr) :: y_d
174  type(c_ptr) :: z_d
175  real(kind=rp), intent(in), optional :: t
176  integer, intent(in), optional :: tstep
177  end subroutine bc_apply_vector_dev
178  end interface
179 
180  interface bc_list_apply
182  end interface bc_list_apply
183 
184  public :: bc_list_init, bc_list_free, bc_list_add, &
186 
187 contains
188 
191  subroutine bc_init(this, coef)
192  class(bc_t), intent(inout) :: this
193  type(coef_t), target, intent(in) :: coef
194 
195  call bc_free(this)
196 
197  this%dof => coef%dof
198  this%coef => coef
199  this%Xh => this%dof%Xh
200  this%msh => this%dof%msh
201 
202  call this%marked_facet%init()
203 
204  end subroutine bc_init
205 
207  subroutine bc_free(this)
208  class(bc_t), intent(inout) :: this
209 
210  call this%marked_facet%free()
211 
212  nullify(this%Xh)
213  nullify(this%msh)
214  nullify(this%dof)
215 
216  if (allocated(this%msk)) then
217  deallocate(this%msk)
218  end if
219 
220  if (allocated(this%facet)) then
221  deallocate(this%facet)
222  end if
223 
224  if (c_associated(this%msk_d)) then
225  call device_free(this%msk_d)
226  this%msk_d = c_null_ptr
227  end if
228 
229  if (c_associated(this%facet_d)) then
230  call device_free(this%facet_d)
231  this%facet_d = c_null_ptr
232  end if
233 
234  end subroutine bc_free
235 
239  subroutine bc_mark_facet(this, facet, el)
240  class(bc_t), intent(inout) :: this
241  integer, intent(in) :: facet
242  integer, intent(in) :: el
243  type(tuple_i4_t) :: t
244 
245  t%x = (/facet, el/)
246  call this%marked_facet%push(t)
247 
248  end subroutine bc_mark_facet
249 
252  subroutine bc_mark_facets(this, facet_list)
253  class(bc_t), intent(inout) :: this
254  type(stack_i4t2_t), intent(inout) :: facet_list
255  type(tuple_i4_t), pointer :: fp(:)
256  integer :: i
257 
258  fp => facet_list%array()
259  do i = 1, facet_list%size()
260  call this%marked_facet%push(fp(i))
261  end do
262 
263  end subroutine bc_mark_facets
264 
267  subroutine bc_mark_zone(this, bc_zone)
268  class(bc_t), intent(inout) :: this
269  class(facet_zone_t), intent(inout) :: bc_zone
270  integer :: i
271  do i = 1, bc_zone%size
272  call this%marked_facet%push(bc_zone%facet_el(i))
273  end do
274  end subroutine bc_mark_zone
275 
282  subroutine bc_mark_zones_from_list(this, bc_zones, bc_key, bc_labels)
283  class(bc_t), intent(inout) :: this
284  class(facet_zone_t), intent(inout) :: bc_zones(:)
285  character(len=*) :: bc_key
286  character(len=100), allocatable :: split_key(:)
287  character(len=NEKO_MSH_MAX_ZLBL_LEN) :: bc_labels(NEKO_MSH_MAX_ZLBLS)
288  integer :: i, j, k, l, msh_bc_type
289 
290  msh_bc_type = 0
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
293  msh_bc_type = 1
294  else if(trim(bc_key) .eq. 'd_pres') then
295  msh_bc_type = 1
296  else if(trim(bc_key) .eq. 'w') then
297  msh_bc_type = 2
298  else if(trim(bc_key) .eq. 'v') then
299  msh_bc_type = 2
300  else if(trim(bc_key) .eq. 'd_vel_u') then
301  msh_bc_type = 2
302  else if(trim(bc_key) .eq. 'd_vel_v') then
303  msh_bc_type = 2
304  else if(trim(bc_key) .eq. 'd_vel_w') then
305  msh_bc_type = 2
306  else if(trim(bc_key) .eq. 'sym') then
307  msh_bc_type = 2
308  end if
309 
310  do i = 1, neko_msh_max_zlbls
311  !Check if several bcs are defined for this zone
312  !bcs are seperated by /, but we could use something else
313  if (index(trim(bc_labels(i)), '/') .eq. 0) then
314  if (trim(bc_key) .eq. trim(bc_labels(i))) then
315  call bc_mark_zone(this, bc_zones(i))
316  ! Loop across all faces in the mesh
317  do j = 1,this%msh%nelv
318  do k = 1, 2 * this%msh%gdim
319  if (this%msh%facet_type(k,j) .eq. -i) then
320  this%msh%facet_type(k,j) = msh_bc_type
321  end if
322  end do
323  end do
324  end if
325  else
326  split_key = split_string(trim(bc_labels(i)),'/')
327  do l = 1, size(split_key)
328  if (trim(split_key(l)) .eq. trim(bc_key)) then
329  call bc_mark_zone(this, bc_zones(i))
330  ! Loop across all faces in the mesh
331  do j = 1,this%msh%nelv
332  do k = 1, 2 * this%msh%gdim
333  if (this%msh%facet_type(k,j) .eq. -i) then
334  this%msh%facet_type(k,j) = msh_bc_type
335  end if
336  end do
337  end do
338  end if
339  end do
340  end if
341  end do
342  end subroutine bc_mark_zones_from_list
343 
344 
348  subroutine bc_finalize(this)
349  class(bc_t), target, intent(inout) :: this
350  type(tuple_i4_t), pointer :: bfp(:)
351  type(tuple_i4_t) :: bc_facet
352  integer :: facet_size, facet, el
353  integer :: i, j, k, l, msk_c
354  integer :: lx, ly, lz, n
355 
356  lx = this%Xh%lx
357  ly = this%Xh%ly
358  lz = this%Xh%lz
359 
361 
362  ! Note we assume that lx = ly = lz
363  facet_size = lx**2
364  allocate(this%msk(0:facet_size * this%marked_facet%size()))
365  allocate(this%facet(0:facet_size * this%marked_facet%size()))
366 
367  msk_c = 0
368  bfp => this%marked_facet%array()
369 
370  ! Loop through each (facet, element) id tuple
371  ! Then loop over all the nodes of the face and compute their linear index
372  ! This index goes into this%msk, whereas the corresponding face id goes into
373  ! this%facet
374  do i = 1, this%marked_facet%size()
375  bc_facet = bfp(i)
376  facet = bc_facet%x(1)
377  el = bc_facet%x(2)
378  select case (facet)
379  case (1)
380  do l = 1, lz
381  do k = 1, ly
382  msk_c = msk_c + 1
383  this%msk(msk_c) = linear_index(1,k,l,el,lx,ly,lz)
384  this%facet(msk_c) = 1
385  end do
386  end do
387  case (2)
388  do l = 1, lz
389  do k = 1, ly
390  msk_c = msk_c + 1
391  this%msk(msk_c) = linear_index(lx,k,l,el,lx,ly,lz)
392  this%facet(msk_c) = 2
393  end do
394  end do
395  case(3)
396  do l = 1, lz
397  do j = 1, lx
398  msk_c = msk_c + 1
399  this%msk(msk_c) = linear_index(j,1,l,el,lx,ly,lz)
400  this%facet(msk_c) = 3
401  end do
402  end do
403  case(4)
404  do l = 1, lz
405  do j = 1, lx
406  msk_c = msk_c + 1
407  this%msk(msk_c) = linear_index(j,ly,l,el,lx,ly,lz)
408  this%facet(msk_c) = 4
409  end do
410  end do
411  case(5)
412  do k = 1, ly
413  do j = 1, lx
414  msk_c = msk_c + 1
415  this%msk(msk_c) = linear_index(j,k,1,el,lx,ly,lz)
416  this%facet(msk_c) = 5
417  end do
418  end do
419  case(6)
420  do k = 1, ly
421  do j = 1, lx
422  msk_c = msk_c + 1
423  this%msk(msk_c) = linear_index(j,k,lz,el,lx,ly,lz)
424  this%facet(msk_c) = 6
425  end do
426  end do
427  end select
428  end do
429 
430  this%msk(0) = msk_c
431  this%facet(0) = msk_c
432 
433  if (neko_bcknd_device .eq. 1) then
434  n = facet_size * this%marked_facet%size() + 1
435  call device_map(this%msk, this%msk_d, n)
436  call device_map(this%facet, this%facet_d, n)
437 
438  call device_memcpy(this%msk, this%msk_d, n, &
439  host_to_device, sync=.false.)
440  call device_memcpy(this%facet, this%facet_d, n, &
441  host_to_device, sync=.false.)
442  end if
443 
444  end subroutine bc_finalize
445 
448  subroutine bc_list_init(bclst, size)
449  type(bc_list_t), intent(inout), target :: bclst
450  integer, optional :: size
451  integer :: n, i
452 
453  call bc_list_free(bclst)
454 
455  if (present(size)) then
456  n = size
457  else
458  n = 1
459  end if
460 
461  allocate(bclst%bc(n))
462 
463  do i = 1, n
464  bclst%bc(i)%bcp => null()
465  end do
466 
467  bclst%n = 0
468  bclst%size = n
469 
470  end subroutine bc_list_init
471 
475  subroutine bc_list_free(bclst)
476  type(bc_list_t), intent(inout) :: bclst
477 
478  if (allocated(bclst%bc)) then
479  deallocate(bclst%bc)
480  end if
481 
482  bclst%n = 0
483  bclst%size = 0
484 
485  end subroutine bc_list_free
486 
489  subroutine bc_list_add(bclst, bc)
490  type(bc_list_t), intent(inout) :: bclst
491  class(bc_t), intent(inout), target :: bc
492  type(bcp_t), allocatable :: tmp(:)
493 
495  if(bc%marked_facet%size() .eq. 0) return
496 
497  if (bclst%n .ge. bclst%size) then
498  bclst%size = bclst%size * 2
499  allocate(tmp(bclst%size))
500  tmp(1:bclst%n) = bclst%bc
501  call move_alloc(tmp, bclst%bc)
502  end if
503 
504  bclst%n = bclst%n + 1
505  bclst%bc(bclst%n)%bcp => bc
506 
507  end subroutine bc_list_add
508 
514  subroutine bc_list_apply_scalar(bclst, x, n, t, tstep)
515  type(bc_list_t), intent(inout) :: bclst
516  integer, intent(in) :: n
517  real(kind=rp), intent(inout), dimension(n) :: x
518  real(kind=rp), intent(in), optional :: t
519  integer, intent(in), optional :: tstep
520  type(c_ptr) :: x_d
521  integer :: i
522 
523  if (neko_bcknd_device .eq. 1) then
524  x_d = device_get_ptr(x)
525  if (present(t) .and. present(tstep)) then
526  do i = 1, bclst%n
527  call bclst%bc(i)%bcp%apply_scalar_dev(x_d, t=t, tstep=tstep)
528  end do
529  else if (present(t)) then
530  do i = 1, bclst%n
531  call bclst%bc(i)%bcp%apply_scalar_dev(x_d, t=t)
532  end do
533  else if (present(tstep)) then
534  do i = 1, bclst%n
535  call bclst%bc(i)%bcp%apply_scalar_dev(x_d, tstep=tstep)
536  end do
537  else
538  do i = 1, bclst%n
539  call bclst%bc(i)%bcp%apply_scalar_dev(x_d)
540  end do
541  end if
542  else
543  if (present(t) .and. present(tstep)) then
544  do i = 1, bclst%n
545  call bclst%bc(i)%bcp%apply_scalar(x, n, t, tstep)
546  end do
547  else if (present(t)) then
548  do i = 1, bclst%n
549  call bclst%bc(i)%bcp%apply_scalar(x, n, t=t)
550  end do
551  else if (present(tstep)) then
552  do i = 1, bclst%n
553  call bclst%bc(i)%bcp%apply_scalar(x, n, tstep=tstep)
554  end do
555  else
556  do i = 1, bclst%n
557  call bclst%bc(i)%bcp%apply_scalar(x, n)
558  end do
559  end if
560  end if
561 
562  end subroutine bc_list_apply_scalar
563 
571  subroutine bc_list_apply_vector(bclst, x, y, z, n, t, tstep)
572  type(bc_list_t), intent(inout) :: bclst
573  integer, intent(in) :: n
574  real(kind=rp), intent(inout), dimension(n) :: x
575  real(kind=rp), intent(inout), dimension(n) :: y
576  real(kind=rp), intent(inout), dimension(n) :: z
577  real(kind=rp), intent(in), optional :: t
578  integer, intent(in), optional :: tstep
579  type(c_ptr) :: x_d
580  type(c_ptr) :: y_d
581  type(c_ptr) :: z_d
582  integer :: i
583 
584  if (neko_bcknd_device .eq. 1) then
585  x_d = device_get_ptr(x)
586  y_d = device_get_ptr(y)
587  z_d = device_get_ptr(z)
588  if (present(t) .and. present(tstep)) then
589  do i = 1, bclst%n
590  call bclst%bc(i)%bcp%apply_vector_dev(x_d, y_d, z_d, t, tstep)
591  end do
592  else if (present(t)) then
593  do i = 1, bclst%n
594  call bclst%bc(i)%bcp%apply_vector_dev(x_d, y_d, z_d, t=t)
595  end do
596  else if (present(tstep)) then
597  do i = 1, bclst%n
598  call bclst%bc(i)%bcp%apply_vector_dev(x_d, y_d, z_d, tstep=tstep)
599  end do
600  else
601  do i = 1, bclst%n
602  call bclst%bc(i)%bcp%apply_vector_dev(x_d, y_d, z_d)
603  end do
604  end if
605  else
606  if (present(t) .and. present(tstep)) then
607  do i = 1, bclst%n
608  call bclst%bc(i)%bcp%apply_vector(x, y, z, n, t, tstep)
609  end do
610  else if (present(t)) then
611  do i = 1, bclst%n
612  call bclst%bc(i)%bcp%apply_vector(x, y, z, n, t=t)
613  end do
614  else if (present(tstep)) then
615  do i = 1, bclst%n
616  call bclst%bc(i)%bcp%apply_vector(x, y, z, n, tstep=tstep)
617  end do
618  else
619  do i = 1, bclst%n
620  call bclst%bc(i)%bcp%apply_vector(x, y, z, n)
621  end do
622  end if
623  end if
624 
625  end subroutine bc_list_apply_vector
626 
627 
628 end module bc
Apply the boundary condition to a scalar field on the device.
Definition: bc.f90:151
Apply the boundary condition to a scalar field.
Definition: bc.f90:116
Apply the boundary condition to a vector field on the device.
Definition: bc.f90:167
Apply the boundary condition to a vector field.
Definition: bc.f90:135
Return the device pointer for an associated Fortran array.
Definition: device.F90:81
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
Copy data between host and device (or device and device)
Definition: device.F90:51
Defines a boundary condition.
Definition: bc.f90:34
subroutine bc_mark_zone(this, bc_zone)
Mark all facets from a zone.
Definition: bc.f90:268
subroutine, public bc_list_add(bclst, bc)
Add a condition to a list of boundary conditions.
Definition: bc.f90:490
subroutine, public bc_list_apply_vector(bclst, x, y, z, n, t, tstep)
Apply a list of boundary conditions to a vector field.
Definition: bc.f90:572
subroutine bc_finalize(this)
Finalize the construction of the bc by populting the msk and facet arrays.
Definition: bc.f90:349
subroutine, public bc_list_init(bclst, size)
Constructor for a list of boundary conditions.
Definition: bc.f90:449
subroutine bc_free(this)
Destructor.
Definition: bc.f90:208
subroutine, public bc_list_free(bclst)
Destructor for a list of boundary conditions.
Definition: bc.f90:476
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...
Definition: bc.f90:283
subroutine bc_mark_facet(this, facet, el)
Mark facet on element el as part of the boundary condition.
Definition: bc.f90:240
subroutine bc_mark_facets(this, facet_list)
Mark all facets from a (facet, el) tuple list.
Definition: bc.f90:253
subroutine, public bc_list_apply_scalar(bclst, x, n, t, tstep)
Apply a list of boundary conditions to a scalar field.
Definition: bc.f90:515
subroutine bc_init(this, coef)
Constructor.
Definition: bc.f90:192
Coefficients.
Definition: coef.f90:34
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
integer, parameter, public host_to_device
Definition: device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition: device.F90:172
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Defines a zone as a subset of facets in a mesh.
Definition: facet_zone.f90:34
Defines a mesh.
Definition: mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition: mesh.f90:55
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition: mesh.f90:57
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a function space.
Definition: space.f90:34
Implements a dynamic stack ADT.
Definition: stack.f90:35
Implements a n-tuple.
Definition: tuple.f90:34
Utilities.
Definition: utils.f90:35
character(len=100) function, dimension(:), allocatable split_string(string, delimiter)
Split a string based on delimiter (tokenizer) OBS: very hacky, this should really be improved,...
Definition: utils.f90:82
pure integer function 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, :)
Definition: utils.f90:123
A list of boundary conditions.
Definition: bc.f90:102
Base type for a boundary condition.
Definition: bc.f90:51
Pointer to boundary condtiion.
Definition: bc.f90:97
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
The function space for the SEM solution fields.
Definition: space.f90:62
Integer 2-tuple based stack.
Definition: stack.f90:84
Integer based 2-tuple.
Definition: tuple.f90:51