Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
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!
34module bc
36 use num_types, only : rp
39 use iso_c_binding, only: c_associated
40 use dofmap, only : dofmap_t
41 use coefs, only : coef_t
42 use space, only : space_t
44 use facet_zone, only : facet_zone_t
45 use stack, only : stack_i4t2_t
46 use tuple, only : tuple_i4_t
47 use field, only : field_t
49 use math, only : relcmp
51 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
52 use json_module, only : json_file
53 implicit none
54 private
55
57 type, public, abstract :: bc_t
59 integer, allocatable :: msk(:)
61 integer, allocatable :: facet(:)
63 type(dofmap_t), pointer :: dof => null()
65 type(coef_t), pointer :: coef => null()
67 type(mesh_t), pointer :: msh => null()
69 type(space_t), pointer :: xh => null()
71 type(stack_i4t2_t) :: marked_facet
73 type(c_ptr) :: msk_d = c_null_ptr
75 type(c_ptr) :: facet_d = c_null_ptr
80 logical :: strong = .true.
81 contains
83 procedure, pass(this) :: init_base => bc_init_base
85 procedure, pass(this) :: free_base => bc_free_base
87 procedure, pass(this) :: mark_facet => bc_mark_facet
89 procedure, pass(this) :: mark_facets => bc_mark_facets
91 procedure, pass(this) :: mark_zone => bc_mark_zone
94 procedure, pass(this) :: finalize_base => bc_finalize_base
95
98 procedure, pass(this) :: apply_scalar_generic => bc_apply_scalar_generic
101 procedure, pass(this) :: apply_vector_generic => bc_apply_vector_generic
103 procedure, pass(this) :: debug_mask_ => bc_debug_mask
105 procedure(bc_apply_scalar), pass(this), deferred :: apply_scalar
107 procedure(bc_apply_vector), pass(this), deferred :: apply_vector
109 procedure(bc_apply_scalar_dev), pass(this), deferred :: apply_scalar_dev
111 procedure(bc_apply_vector_dev), pass(this), deferred :: apply_vector_dev
113 procedure(bc_destructor), pass(this), deferred :: free
115 procedure(bc_constructor), pass(this), deferred :: init
117 procedure(bc_finalize), pass(this), deferred :: finalize
118 end type bc_t
119
121 type, public :: bc_ptr_t
122 class(bc_t), pointer :: ptr => null()
123 end type bc_ptr_t
124
125 ! Helper type to have an array of polymorphic bc_t objects.
126 type, public :: bc_alloc_t
127 class(bc_t), allocatable :: obj
128 end type bc_alloc_t
129
130
131 abstract interface
132
133 subroutine bc_constructor(this, coef, json)
134 import :: bc_t, coef_t, json_file
135 class(bc_t), intent(inout), target :: this
136 type(coef_t), intent(in) :: coef
137 type(json_file), intent(inout) :: json
138 end subroutine bc_constructor
139 end interface
140
141 abstract interface
142
143 subroutine bc_destructor(this)
144 import :: bc_t
145 class(bc_t), intent(inout), target :: this
146 end subroutine bc_destructor
147 end interface
148
149 abstract interface
150
151 subroutine bc_finalize(this)
152 import :: bc_t
153 class(bc_t), intent(inout), target :: this
154 end subroutine bc_finalize
155 end interface
156
157 abstract interface
158
164 subroutine bc_apply_scalar(this, x, n, t, tstep, strong)
165 import :: bc_t
166 import :: rp
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
173 end subroutine bc_apply_scalar
174 end interface
175
176 abstract interface
177
185 subroutine bc_apply_vector(this, x, y, z, n, t, tstep, strong)
186 import :: bc_t
187 import :: rp
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
196 end subroutine bc_apply_vector
197 end interface
198
199 abstract interface
200
205 subroutine bc_apply_scalar_dev(this, x_d, t, tstep, strong)
206 import :: c_ptr
207 import :: bc_t
208 import :: rp
209 class(bc_t), intent(inout), target :: this
210 type(c_ptr) :: x_d
211 real(kind=rp), intent(in), optional :: t
212 integer, intent(in), optional :: tstep
213 logical, intent(in), optional :: strong
214 end subroutine bc_apply_scalar_dev
215 end interface
216
217 abstract interface
218
225 subroutine bc_apply_vector_dev(this, x_d, y_d, z_d, t, tstep, strong)
226 import :: c_ptr
227 import :: bc_t
228 import :: rp
229 class(bc_t), intent(inout), target :: this
230 type(c_ptr) :: x_d
231 type(c_ptr) :: y_d
232 type(c_ptr) :: z_d
233 real(kind=rp), intent(in), optional :: t
234 integer, intent(in), optional :: tstep
235 logical, intent(in), optional :: strong
236 end subroutine bc_apply_vector_dev
237 end interface
238
239contains
240
243 subroutine bc_init_base(this, coef)
244 class(bc_t), intent(inout) :: this
245 type(coef_t), target, intent(in) :: coef
246
247 call this%free_base
248
249 this%dof => coef%dof
250 this%coef => coef
251 this%Xh => this%dof%Xh
252 this%msh => this%dof%msh
253
254 call this%marked_facet%init()
255
256 end subroutine bc_init_base
257
259 subroutine bc_free_base(this)
260 class(bc_t), intent(inout) :: this
261
262 call this%marked_facet%free()
263
264 nullify(this%Xh)
265 nullify(this%msh)
266 nullify(this%dof)
267 nullify(this%coef)
268
269 if (allocated(this%msk)) then
270 deallocate(this%msk)
271 end if
272
273 if (allocated(this%facet)) then
274 deallocate(this%facet)
275 end if
276
277 if (c_associated(this%msk_d)) then
278 call device_free(this%msk_d)
279 this%msk_d = c_null_ptr
280 end if
281
282 if (c_associated(this%facet_d)) then
283 call device_free(this%facet_d)
284 this%facet_d = c_null_ptr
285 end if
286
287 end subroutine bc_free_base
288
297 subroutine bc_apply_vector_generic(this, x, y, z, n, t, tstep)
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
305 type(c_ptr) :: x_d
306 type(c_ptr) :: y_d
307 type(c_ptr) :: z_d
308 integer :: i
309
310
311 if (neko_bcknd_device .eq. 1) then
312 x_d = device_get_ptr(x)
313 y_d = device_get_ptr(y)
314 z_d = device_get_ptr(z)
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)
321 else
322 call this%apply_vector_dev(x_d, y_d, z_d)
323 end if
324 else
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)
331 else
332 call this%apply_vector(x, y, z, n)
333 end if
334 end if
335
336 end subroutine bc_apply_vector_generic
337
346 subroutine bc_apply_scalar_generic(this, x, n, t, tstep)
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
352 type(c_ptr) :: x_d
353 integer :: i
354
355
356 if (neko_bcknd_device .eq. 1) then
357 x_d = device_get_ptr(x)
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)
364 else
365 call this%apply_scalar_dev(x_d)
366 end if
367 else
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)
374 else
375 call this%apply_scalar(x, n)
376 end if
377 end if
378
379 end subroutine bc_apply_scalar_generic
380
384 subroutine bc_mark_facet(this, facet, el)
385 class(bc_t), intent(inout) :: this
386 integer, intent(in) :: facet
387 integer, intent(in) :: el
388 type(tuple_i4_t) :: t
389
390 t%x = [facet, el]
391 call this%marked_facet%push(t)
392
393 end subroutine bc_mark_facet
394
397 subroutine bc_mark_facets(this, facet_list)
398 class(bc_t), intent(inout) :: this
399 type(stack_i4t2_t), intent(inout) :: facet_list
400 type(tuple_i4_t), pointer :: fp(:)
401 integer :: i
402
403 fp => facet_list%array()
404 do i = 1, facet_list%size()
405 call this%marked_facet%push(fp(i))
406 end do
407
408 end subroutine bc_mark_facets
409
412 subroutine bc_mark_zone(this, bc_zone)
413 class(bc_t), intent(inout) :: this
414 class(facet_zone_t), intent(in) :: bc_zone
415 integer :: i
416 do i = 1, bc_zone%size
417 call this%marked_facet%push(bc_zone%facet_el(i))
418 end do
419 end subroutine bc_mark_zone
420
428 subroutine bc_finalize_base(this, only_facets)
429 class(bc_t), target, intent(inout) :: this
430 logical, optional, intent(in) :: only_facets
431 type(tuple_i4_t), pointer :: bfp(:)
432 type(tuple_i4_t) :: bc_facet
433 type(field_t) :: test_field
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
438 lx = this%Xh%lx
439 ly = this%Xh%ly
440 lz = this%Xh%lz
441 if ( present(only_facets)) then
442 only_facet = only_facets
443 else
444 only_facet = .false.
445 end if
447
448 ! Note we assume that lx = ly = lz
449 facet_size = lx**2
450 allocate(this%msk(0:facet_size * this%marked_facet%size()))
451 allocate(this%facet(0:facet_size * this%marked_facet%size()))
452
453 msk_c = 0
454 bfp => this%marked_facet%array()
455
456 ! Loop through each (facet, element) id tuple
457 ! Then loop over all the nodes of the face and compute their linear index
458 ! This index goes into this%msk, whereas the corresponding face id goes into
459 ! this%facet
460 do i = 1, this%marked_facet%size()
461 bc_facet = bfp(i)
462 facet = bc_facet%x(1)
463 el = bc_facet%x(2)
464 select case (facet)
465 case (1)
466 do l = 1, lz
467 do k = 1, ly
468 msk_c = msk_c + 1
469 this%msk(msk_c) = linear_index(1, k, l, el, lx, ly, lz)
470 this%facet(msk_c) = 1
471 end do
472 end do
473 case (2)
474 do l = 1, lz
475 do k = 1, ly
476 msk_c = msk_c + 1
477 this%msk(msk_c) = linear_index(lx, k, l, el, lx, ly, lz)
478 this%facet(msk_c) = 2
479 end do
480 end do
481 case (3)
482 do l = 1, lz
483 do j = 1, lx
484 msk_c = msk_c + 1
485 this%msk(msk_c) = linear_index(j, 1, l, el, lx, ly, lz)
486 this%facet(msk_c) = 3
487 end do
488 end do
489 case (4)
490 do l = 1, lz
491 do j = 1, lx
492 msk_c = msk_c + 1
493 this%msk(msk_c) = linear_index(j, ly, l, el, lx, ly, lz)
494 this%facet(msk_c) = 4
495 end do
496 end do
497 case (5)
498 do k = 1, ly
499 do j = 1, lx
500 msk_c = msk_c + 1
501 this%msk(msk_c) = linear_index(j, k, 1, el, lx, ly, lz)
502 this%facet(msk_c) = 5
503 end do
504 end do
505 case (6)
506 do k = 1, ly
507 do j = 1, lx
508 msk_c = msk_c + 1
509 this%msk(msk_c) = linear_index(j, k, lz, el, lx, ly, lz)
510 this%facet(msk_c) = 6
511 end do
512 end do
513 end select
514 end do
515 if ( .not. only_facet) then
516 !Makes check for points not on facet that should have bc applied
517 call test_field%init(this%dof)
518
519 n = test_field%size()
520 test_field%x = 0.0_rp
521 !Apply this bc once
522 do i = 1, msk_c
523 test_field%x(this%msk(i),1,1,1) = 1.0
524 end do
525 if (neko_bcknd_device .eq. 1) then
526 call device_memcpy(test_field%x, test_field%x_d, n, &
527 host_to_device, sync=.true.)
528 end if
529 !Check if some point that was not zeroed was zeroed on another element
530 call this%coef%gs_h%op(test_field,gs_op_add)
531 if (neko_bcknd_device .eq. 1) then
532 call device_memcpy(test_field%x, test_field%x_d, n, &
533 device_to_host, sync=.true.)
534 end if
535 msk_c = 0
536 do i = 1, this%dof%size()
537 if (test_field%x(i,1,1,1) .gt. 0.5) then
538 msk_c = msk_c + 1
539 end if
540 end do
541 !Allocate new mask
542 deallocate(this%msk)
543 allocate(this%msk(0:msk_c))
544 j = 1
545 do i = 1, this%dof%size()
546 if (test_field%x(i,1,1,1) .gt. 0.5) then
547 this%msk(j) = i
548 j = j + 1
549 end if
550 end do
551
552 call test_field%free()
553 end if
554
555 this%msk(0) = msk_c
556 this%facet(0) = msk_c
557
558 if (neko_bcknd_device .eq. 1) then
559 n = msk_c + 1
560 call device_map(this%msk, this%msk_d, n)
561 call device_map(this%facet, this%facet_d, n)
562
563 call device_memcpy(this%msk, this%msk_d, n, &
564 host_to_device, sync = .false.)
565 call device_memcpy(this%facet, this%facet_d, n, &
566 host_to_device, sync = .true.)
567 end if
568
569 end subroutine bc_finalize_base
570
574 subroutine bc_debug_mask(this, file_name)
575 use field, only : field_t
576 use file, only : file_t
577
578 class(bc_t), intent(inout) :: this
579 character(len=*), intent(in) :: file_name
580 type(field_t) :: bdry_field
581 integer:: i, m, k
582 type(file_t) :: dump_file
583
584 call bdry_field%init(this%coef%dof, 'bdry')
585 m = this%msk(0)
586 do i = 1, m
587 k = this%msk(i)
588 bdry_field%x(k,1,1,1) = 1.0_rp
589 end do
590 dump_file = file_t(file_name)
591 call dump_file%write(bdry_field)
592
593 end subroutine bc_debug_mask
594end module bc
Apply the boundary condition to a scalar field on the device.
Definition bc.f90:205
Apply the boundary condition to a scalar field.
Definition bc.f90:164
Apply the boundary condition to a vector field on the device.
Definition bc.f90:225
Apply the boundary condition to a vector field.
Definition bc.f90:185
Constructor.
Definition bc.f90:133
Destructor.
Definition bc.f90:143
Finalize by building the mask and facet arrays.
Definition bc.f90:151
Return the device pointer for an associated Fortran array.
Definition device.F90:95
Map a Fortran array to a device (allocate and associate)
Definition device.F90:71
Copy data between host and device (or device and device)
Definition device.F90:65
Defines a boundary condition.
Definition bc.f90:34
subroutine bc_mark_zone(this, bc_zone)
Mark all facets from a zone.
Definition bc.f90:413
subroutine bc_free_base(this)
Destructor for the base type, bc_t.
Definition bc.f90:260
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.
Definition bc.f90:347
subroutine bc_init_base(this, coef)
Constructor.
Definition bc.f90:244
subroutine bc_finalize_base(this, only_facets)
Finalize the construction of the bc by populting the msk and facet arrays.
Definition bc.f90:429
subroutine bc_debug_mask(this, file_name)
Write a field showing the mask of the bc.
Definition bc.f90:575
subroutine bc_mark_facet(this, facet, el)
Mark facet on element el as part of the boundary condition.
Definition bc.f90:385
subroutine bc_mark_facets(this, facet_list)
Mark all facets from a (facet, el) tuple list.
Definition bc.f90:398
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.
Definition bc.f90:298
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:46
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:200
integer, parameter, public device_to_host
Definition device.F90:46
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a zone as a subset of facets in a mesh.
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Gather-scatter.
Definition math.f90:60
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:56
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:58
Build configurations.
integer, parameter neko_bcknd_device
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, public split_string(string, delimiter)
Split a string based on delimiter (tokenizer) OBS: very hacky, this should really be improved,...
Definition utils.f90:136
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, :)
Definition utils.f90:174
Pointer to a `bc_t`.
Definition bc.f90:121
Base type for a boundary condition.
Definition bc.f90:57
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.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