Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
bc.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2025, 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
37 use device, only : host_to_device, device_memcpy, &
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
48 use gs_ops, only : gs_op_add
49 use math, only : relcmp
51 use logger, only : neko_log, log_size
52 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
53 use json_module, only : json_file
54 use time_state, only : time_state_t
55 use field, only : field_t
56 use file, only : file_t
57
58 implicit none
59 private
60
62 type, public, abstract :: bc_t
64 integer, allocatable :: msk(:)
66 integer, allocatable :: facet(:)
68 type(dofmap_t), pointer :: dof => null()
70 type(coef_t), pointer :: coef => null()
72 type(mesh_t), pointer :: msh => null()
74 type(space_t), pointer :: xh => null()
76 type(stack_i4t2_t) :: marked_facet
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.
89 !!> Name of the bc
90 character(len=:), allocatable :: name
91 !!> Zone indices where the bc is applied
92 integer, allocatable :: zone_indices(:)
93 contains
95 procedure, pass(this) :: init_base => bc_init_base
97 procedure, pass(this) :: free_base => bc_free_base
99 procedure, pass(this) :: mark_facet => bc_mark_facet
101 procedure, pass(this) :: mark_facets => bc_mark_facets
103 procedure, pass(this) :: mark_zone => bc_mark_zone
106 procedure, pass(this) :: finalize_base => bc_finalize_base
107
110 procedure, pass(this) :: apply_scalar_generic => bc_apply_scalar_generic
113 procedure, pass(this) :: apply_vector_generic => bc_apply_vector_generic
115 procedure, pass(this) :: debug_mask_ => bc_debug_mask
117 procedure(bc_apply_scalar), pass(this), deferred :: apply_scalar
119 procedure(bc_apply_vector), pass(this), deferred :: apply_vector
121 procedure(bc_apply_scalar_dev), pass(this), deferred :: apply_scalar_dev
123 procedure(bc_apply_vector_dev), pass(this), deferred :: apply_vector_dev
125 procedure(bc_destructor), pass(this), deferred :: free
127 procedure(bc_constructor), pass(this), deferred :: init
129 procedure(bc_finalize), pass(this), deferred :: finalize
130 end type bc_t
131
133 type, public :: bc_ptr_t
134 class(bc_t), pointer :: ptr => null()
135 end type bc_ptr_t
136
137 ! Helper type to have an array of polymorphic bc_t objects.
138 type, public :: bc_alloc_t
139 class(bc_t), allocatable :: obj
140 end type bc_alloc_t
141
142
143 abstract interface
144
145 subroutine bc_constructor(this, coef, json)
146 import :: bc_t, coef_t, json_file
147 class(bc_t), intent(inout), target :: this
148 type(coef_t), target, intent(in) :: coef
149 type(json_file), intent(inout) :: json
150 end subroutine bc_constructor
151 end interface
152
153 abstract interface
154
155 subroutine bc_destructor(this)
156 import :: bc_t
157 class(bc_t), intent(inout), target :: this
158 end subroutine bc_destructor
159 end interface
160
161 abstract interface
162
163 subroutine bc_finalize(this, only_facets)
164 import :: bc_t
165 class(bc_t), intent(inout), target :: this
166 logical, optional, intent(in) :: only_facets
167 end subroutine bc_finalize
168 end interface
169
170 abstract interface
171
176 subroutine bc_apply_scalar(this, x, n, time, strong)
177 import :: bc_t, time_state_t
178 import :: rp
179 class(bc_t), intent(inout) :: this
180 integer, intent(in) :: n
181 real(kind=rp), intent(inout), dimension(n) :: x
182 type(time_state_t), intent(in), optional :: time
183 logical, intent(in), optional :: strong
184 end subroutine bc_apply_scalar
185 end interface
186
187 abstract interface
188
196 subroutine bc_apply_vector(this, x, y, z, n, time, strong)
197 import :: bc_t, time_state_t
198 import :: rp
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
204 type(time_state_t), intent(in), optional :: time
205 logical, intent(in), optional :: strong
206 end subroutine bc_apply_vector
207 end interface
208
209 abstract interface
210
215 subroutine bc_apply_scalar_dev(this, x_d, time, strong, strm)
216 import :: c_ptr
217 import :: bc_t, time_state_t
218 import :: rp
219 class(bc_t), intent(inout), target :: this
220 type(c_ptr), intent(inout) :: x_d
221 type(time_state_t), intent(in), optional :: time
222 logical, intent(in), optional :: strong
223 type(c_ptr), intent(inout) :: strm
224 end subroutine bc_apply_scalar_dev
225 end interface
226
227 abstract interface
228
235 subroutine bc_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
236 import :: c_ptr, bc_t, time_state_t
237 import :: rp
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
242 type(time_state_t), intent(in), optional :: time
243 logical, intent(in), optional :: strong
244 type(c_ptr), intent(inout) :: strm
245 end subroutine bc_apply_vector_dev
246 end interface
247
248contains
249
252 subroutine bc_init_base(this, coef)
253 class(bc_t), intent(inout) :: this
254 type(coef_t), target, intent(in) :: coef
255
256 call this%free_base
257
258 this%dof => coef%dof
259 this%coef => coef
260 this%Xh => this%dof%Xh
261 this%msh => this%dof%msh
262
263 call this%marked_facet%init()
264
265 end subroutine bc_init_base
266
268 subroutine bc_free_base(this)
269 class(bc_t), intent(inout) :: this
270
271 call this%marked_facet%free()
272
273 nullify(this%Xh)
274 nullify(this%msh)
275 nullify(this%dof)
276 nullify(this%coef)
277
278 if (allocated(this%msk)) then
279 if (neko_bcknd_device .eq. 1) then
280 call device_unmap(this%msk, this%msk_d)
281 end if
282 deallocate(this%msk)
283 end if
284
285 if (allocated(this%facet)) then
286 if (neko_bcknd_device .eq. 1) then
287 call device_unmap(this%facet, this%facet_d)
288 end if
289 deallocate(this%facet)
290 end if
291
292 if (allocated(this%name)) then
293 deallocate(this%name)
294 end if
295
296 if (allocated(this%zone_indices)) then
297 deallocate(this%zone_indices)
298 end if
299
300 end subroutine bc_free_base
301
310 subroutine bc_apply_vector_generic(this, x, y, z, time, strong, strm)
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
315 type(time_state_t), intent(in), optional :: time
316 logical, intent(in), optional :: strong
317 type(c_ptr), intent(inout), optional :: strm
318 type(c_ptr) :: strm_
319 integer :: n
320 character(len=256) :: msg
321
322 ! Get the size of the fields
323 n = x%size()
324
325 ! Ensure all fields are the same size
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"
329 call neko_error(trim(msg))
330 end if
331
332 if (neko_bcknd_device .eq. 1) then
333
334 if (present(strm)) then
335 strm_ = strm
336 else
337 strm_ = glb_cmd_queue
338 end if
339
340 call this%apply_vector_dev(x%x_d, y%x_d, z%x_d, time = time, &
341 strong = strong, strm = strm_)
342 else
343 call this%apply_vector(x%x, y%x, z%x, n, time = time, strong = strong)
344 end if
345
346 end subroutine bc_apply_vector_generic
347
354 subroutine bc_apply_scalar_generic(this, x, time, strong, strm)
355 class(bc_t), intent(inout) :: this
356 type(field_t), intent(inout) :: x
357 type(time_state_t), intent(in), optional :: time
358 logical, intent(in), optional :: strong
359 type(c_ptr), intent(inout), optional :: strm
360 type(c_ptr) :: strm_
361 integer :: n
362
363 ! Get the size of the field
364 n = x%size()
365
366 if (neko_bcknd_device .eq. 1) then
367
368 if (present(strm)) then
369 strm_ = strm
370 else
371 strm_ = glb_cmd_queue
372 end if
373
374 call this%apply_scalar_dev(x%x_d, time = time, strong = strong, &
375 strm = strm_)
376 else
377 call this%apply_scalar(x%x, n, time = time)
378 end if
379
380 end subroutine bc_apply_scalar_generic
381
385 subroutine bc_mark_facet(this, facet, el)
386 class(bc_t), intent(inout) :: this
387 integer, intent(in) :: facet
388 integer, intent(in) :: el
389 type(tuple_i4_t) :: t
390
391 t%x = [facet, el]
392 call this%marked_facet%push(t)
393
394 end subroutine bc_mark_facet
395
398 subroutine bc_mark_facets(this, facet_list)
399 class(bc_t), intent(inout) :: this
400 type(stack_i4t2_t), intent(inout) :: facet_list
401 type(tuple_i4_t), pointer :: fp(:)
402 integer :: i
403
404 fp => facet_list%array()
405 do i = 1, facet_list%size()
406 call this%marked_facet%push(fp(i))
407 end do
408
409 end subroutine bc_mark_facets
410
413 subroutine bc_mark_zone(this, bc_zone)
414 class(bc_t), intent(inout) :: this
415 class(facet_zone_t), intent(in) :: bc_zone
416 integer :: i
417 do i = 1, bc_zone%size
418 call this%marked_facet%push(bc_zone%facet_el(i))
419 end do
420 end subroutine bc_mark_zone
421
429 subroutine bc_finalize_base(this, only_facets)
430 class(bc_t), target, intent(inout) :: this
431 logical, optional, intent(in) :: only_facets
432 type(tuple_i4_t), pointer :: bfp(:)
433 type(tuple_i4_t) :: bc_facet
434 type(field_t) :: test_field
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
440 lx = this%Xh%lx
441 ly = this%Xh%ly
442 lz = this%Xh%lz
443 if ( present(only_facets)) then
444 only_facet = only_facets
445 else
446 only_facet = .false.
447 end if
449
450 ! Note we assume that lx = ly = lz
451 facet_size = lx**2
452 allocate(this%msk(0:facet_size * this%marked_facet%size()))
453 allocate(this%facet(0:facet_size * this%marked_facet%size()))
454
455 msk_c = 0
456 bfp => this%marked_facet%array()
457
458 ! Loop through each (facet, element) id tuple
459 ! Then loop over all the nodes of the face and compute their linear index
460 ! This index goes into this%msk, whereas the corresponding face id goes into
461 ! this%facet
462 do i = 1, this%marked_facet%size()
463 bc_facet = bfp(i)
464 facet = bc_facet%x(1)
465 el = bc_facet%x(2)
466 select case (facet)
467 case (1)
468 do l = 1, lz
469 do k = 1, ly
470 msk_c = msk_c + 1
471 this%msk(msk_c) = linear_index(1, k, l, el, lx, ly, lz)
472 this%facet(msk_c) = 1
473 end do
474 end do
475 case (2)
476 do l = 1, lz
477 do k = 1, ly
478 msk_c = msk_c + 1
479 this%msk(msk_c) = linear_index(lx, k, l, el, lx, ly, lz)
480 this%facet(msk_c) = 2
481 end do
482 end do
483 case (3)
484 do l = 1, lz
485 do j = 1, lx
486 msk_c = msk_c + 1
487 this%msk(msk_c) = linear_index(j, 1, l, el, lx, ly, lz)
488 this%facet(msk_c) = 3
489 end do
490 end do
491 case (4)
492 do l = 1, lz
493 do j = 1, lx
494 msk_c = msk_c + 1
495 this%msk(msk_c) = linear_index(j, ly, l, el, lx, ly, lz)
496 this%facet(msk_c) = 4
497 end do
498 end do
499 case (5)
500 do k = 1, ly
501 do j = 1, lx
502 msk_c = msk_c + 1
503 this%msk(msk_c) = linear_index(j, k, 1, el, lx, ly, lz)
504 this%facet(msk_c) = 5
505 end do
506 end do
507 case (6)
508 do k = 1, ly
509 do j = 1, lx
510 msk_c = msk_c + 1
511 this%msk(msk_c) = linear_index(j, k, lz, el, lx, ly, lz)
512 this%facet(msk_c) = 6
513 end do
514 end do
515 end select
516 end do
517 this%facet(0) = msk_c
518 if (neko_bcknd_device .eq. 1) then
519 !Observe the facet_mask is junk if only_facet is false
520 n = msk_c + 1
521 call device_map(this%facet, this%facet_d, n)
522 call device_memcpy(this%facet, this%facet_d, n, &
523 host_to_device, sync = .true.)
524 end if
525 if ( .not. only_facet) then
526 !Makes check for points not on facet that should have bc applied
527 call test_field%init(this%dof)
528
529 n = test_field%size()
530 test_field%x = 0.0_rp
531 !Apply this bc once
532 do i = 1, msk_c
533 test_field%x(this%msk(i),1,1,1) = 1.0
534 end do
535 if (neko_bcknd_device .eq. 1) then
536 call device_memcpy(test_field%x, test_field%x_d, n, &
537 host_to_device, sync = .true.)
538 end if
539 !Check if some point that was not zeroed was zeroed on another element
540 call this%coef%gs_h%op(test_field, gs_op_add)
541 if (neko_bcknd_device .eq. 1) then
542 call device_memcpy(test_field%x, test_field%x_d, n, &
543 device_to_host, sync = .true.)
544 end if
545 msk_c = 0
546 do i = 1, this%dof%size()
547 if (test_field%x(i,1,1,1) .gt. 0.5) then
548 msk_c = msk_c + 1
549 end if
550 end do
551 !Allocate new mask
552 deallocate(this%msk)
553 allocate(this%msk(0:msk_c))
554 j = 1
555 do i = 1, this%dof%size()
556 if (test_field%x(i,1,1,1) .gt. 0.5) then
557 this%msk(j) = i
558 j = j + 1
559 end if
560 end do
561
562 call test_field%free()
563 end if
564
565 this%msk(0) = msk_c
566
567 if (neko_bcknd_device .eq. 1) then
568 n = msk_c + 1
569 call device_map(this%msk, this%msk_d, n)
570 call device_memcpy(this%msk, this%msk_d, n, &
571 host_to_device, sync = .true.)
572 end if
573
574 if (.not. allocated(this%name)) then
575! gives plenty of empty info lines during AMR restart
576! this%name = ""
577 else
578 write(log_buf, '(A,A)') 'BC assigned name : ', trim(this%name)
579 call neko_log%message(log_buf)
580 end if
581
582! causes trouble for AMR restart
583! if (.not. allocated(this%zone_indices)) then
584! allocate(this%zone_indices(1))
585! this%zone_indices(1) = -1
586! end if
587
588 end subroutine bc_finalize_base
589
593 subroutine bc_debug_mask(this, file_name)
594 class(bc_t), intent(inout) :: this
595 character(len=*), intent(in) :: file_name
596 type(field_t) :: bdry_field
597 integer:: i, m, k
598 type(file_t) :: dump_file
599
600 call bdry_field%init(this%coef%dof, 'bdry')
601 m = this%msk(0)
602 do i = 1, m
603 k = this%msk(i)
604 bdry_field%x(k,1,1,1) = 1.0_rp
605 end do
606 call dump_file%init(file_name)
607 call dump_file%write(bdry_field)
608
609 end subroutine bc_debug_mask
610end module bc
Apply the boundary condition to a scalar field on the device.
Definition bc.f90:215
Apply the boundary condition to a scalar field.
Definition bc.f90:176
Apply the boundary condition to a vector field on the device.
Definition bc.f90:235
Apply the boundary condition to a vector field.
Definition bc.f90:196
Constructor.
Definition bc.f90:145
Destructor.
Definition bc.f90:155
Finalize by building the mask and facet arrays.
Definition bc.f90:163
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
Unmap a Fortran array from a device (deassociate and free)
Definition device.F90:83
Defines a boundary condition.
Definition bc.f90:34
subroutine bc_mark_zone(this, bc_zone)
Mark all facets from a zone.
Definition bc.f90:414
subroutine bc_free_base(this)
Destructor for the base type, bc_t.
Definition bc.f90:269
subroutine bc_init_base(this, coef)
Constructor.
Definition bc.f90:253
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.
Definition bc.f90:355
subroutine bc_finalize_base(this, only_facets)
Finalize the construction of the bc by populting the msk and facet arrays.
Definition bc.f90:430
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.
Definition bc.f90:311
subroutine bc_debug_mask(this, file_name)
Write a field showing the mask of the bc.
Definition bc.f90:594
subroutine bc_mark_facet(this, facet, el)
Mark facet on element el as part of the boundary condition.
Definition bc.f90:386
subroutine bc_mark_facets(this, facet_list)
Mark all facets from a (facet, el) tuple list.
Definition bc.f90:399
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
integer, parameter, public device_to_host
Definition device.F90:47
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:51
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
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
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:63
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:65
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:49
Module with things related to the simulation time.
Implements a n-tuple.
Definition tuple.f90:41
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:199
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:237
Pointer to a `bc_t`.
Definition bc.f90:133
Base type for a boundary condition.
Definition bc.f90:62
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:61
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:56
The function space for the SEM solution fields.
Definition space.f90:63
Integer 2-tuple based stack.
Definition stack.f90:98
A struct that contains all info about the time, expand as needed.
Integer based 2-tuple.
Definition tuple.f90:58