Neko 1.99.2
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 deallocate(this%msk)
280 end if
281
282 if (allocated(this%facet)) then
283 deallocate(this%facet)
284 end if
285
286 if (c_associated(this%msk_d)) then
287 call device_free(this%msk_d)
288 this%msk_d = c_null_ptr
289 end if
290
291 if (c_associated(this%facet_d)) then
292 call device_free(this%facet_d)
293 this%facet_d = c_null_ptr
294 end if
295
296 if (allocated(this%name)) then
297 deallocate(this%name)
298 end if
299
300 if (allocated(this%zone_indices)) then
301 deallocate(this%zone_indices)
302 end if
303
304 end subroutine bc_free_base
305
314 subroutine bc_apply_vector_generic(this, x, y, z, time, strong, strm)
315 class(bc_t), intent(inout) :: this
316 type(field_t), intent(inout) :: x
317 type(field_t), intent(inout) :: y
318 type(field_t), intent(inout) :: z
319 type(time_state_t), intent(in), optional :: time
320 logical, intent(in), optional :: strong
321 type(c_ptr), intent(inout), optional :: strm
322 type(c_ptr) :: strm_
323 integer :: n
324 character(len=256) :: msg
325
326 ! Get the size of the fields
327 n = x%size()
328
329 ! Ensure all fields are the same size
330 if (y%size() .ne. n .or. z%size() .ne. n) then
331 msg = "Fields x, y, z must have the same size in " // &
332 "bc_list_apply_vector_field"
333 call neko_error(trim(msg))
334 end if
335
336 if (neko_bcknd_device .eq. 1) then
337
338 if (present(strm)) then
339 strm_ = strm
340 else
341 strm_ = glb_cmd_queue
342 end if
343
344 call this%apply_vector_dev(x%x_d, y%x_d, z%x_d, time = time, &
345 strong = strong, strm = strm_)
346 else
347 call this%apply_vector(x%x, y%x, z%x, n, time = time, strong = strong)
348 end if
349
350 end subroutine bc_apply_vector_generic
351
358 subroutine bc_apply_scalar_generic(this, x, time, strong, strm)
359 class(bc_t), intent(inout) :: this
360 type(field_t), intent(inout) :: x
361 type(time_state_t), intent(in), optional :: time
362 logical, intent(in), optional :: strong
363 type(c_ptr), intent(inout), optional :: strm
364 type(c_ptr) :: strm_
365 integer :: n
366
367 ! Get the size of the field
368 n = x%size()
369
370 if (neko_bcknd_device .eq. 1) then
371
372 if (present(strm)) then
373 strm_ = strm
374 else
375 strm_ = glb_cmd_queue
376 end if
377
378 call this%apply_scalar_dev(x%x_d, time = time, strong = strong, &
379 strm = strm_)
380 else
381 call this%apply_scalar(x%x, n, time = time)
382 end if
383
384 end subroutine bc_apply_scalar_generic
385
389 subroutine bc_mark_facet(this, facet, el)
390 class(bc_t), intent(inout) :: this
391 integer, intent(in) :: facet
392 integer, intent(in) :: el
393 type(tuple_i4_t) :: t
394
395 t%x = [facet, el]
396 call this%marked_facet%push(t)
397
398 end subroutine bc_mark_facet
399
402 subroutine bc_mark_facets(this, facet_list)
403 class(bc_t), intent(inout) :: this
404 type(stack_i4t2_t), intent(inout) :: facet_list
405 type(tuple_i4_t), pointer :: fp(:)
406 integer :: i
407
408 fp => facet_list%array()
409 do i = 1, facet_list%size()
410 call this%marked_facet%push(fp(i))
411 end do
412
413 end subroutine bc_mark_facets
414
417 subroutine bc_mark_zone(this, bc_zone)
418 class(bc_t), intent(inout) :: this
419 class(facet_zone_t), intent(in) :: bc_zone
420 integer :: i
421 do i = 1, bc_zone%size
422 call this%marked_facet%push(bc_zone%facet_el(i))
423 end do
424 end subroutine bc_mark_zone
425
433 subroutine bc_finalize_base(this, only_facets)
434 class(bc_t), target, intent(inout) :: this
435 logical, optional, intent(in) :: only_facets
436 type(tuple_i4_t), pointer :: bfp(:)
437 type(tuple_i4_t) :: bc_facet
438 type(field_t) :: test_field
439 integer :: facet_size, facet, el
440 logical :: only_facet = .false.
441 integer :: i, j, k, l, msk_c
442 integer :: lx, ly, lz, n
443 character(len=LOG_SIZE) :: log_buf
444 lx = this%Xh%lx
445 ly = this%Xh%ly
446 lz = this%Xh%lz
447 if ( present(only_facets)) then
448 only_facet = only_facets
449 else
450 only_facet = .false.
451 end if
453
454 ! Note we assume that lx = ly = lz
455 facet_size = lx**2
456 allocate(this%msk(0:facet_size * this%marked_facet%size()))
457 allocate(this%facet(0:facet_size * this%marked_facet%size()))
458
459 msk_c = 0
460 bfp => this%marked_facet%array()
461
462 ! Loop through each (facet, element) id tuple
463 ! Then loop over all the nodes of the face and compute their linear index
464 ! This index goes into this%msk, whereas the corresponding face id goes into
465 ! this%facet
466 do i = 1, this%marked_facet%size()
467 bc_facet = bfp(i)
468 facet = bc_facet%x(1)
469 el = bc_facet%x(2)
470 select case (facet)
471 case (1)
472 do l = 1, lz
473 do k = 1, ly
474 msk_c = msk_c + 1
475 this%msk(msk_c) = linear_index(1, k, l, el, lx, ly, lz)
476 this%facet(msk_c) = 1
477 end do
478 end do
479 case (2)
480 do l = 1, lz
481 do k = 1, ly
482 msk_c = msk_c + 1
483 this%msk(msk_c) = linear_index(lx, k, l, el, lx, ly, lz)
484 this%facet(msk_c) = 2
485 end do
486 end do
487 case (3)
488 do l = 1, lz
489 do j = 1, lx
490 msk_c = msk_c + 1
491 this%msk(msk_c) = linear_index(j, 1, l, el, lx, ly, lz)
492 this%facet(msk_c) = 3
493 end do
494 end do
495 case (4)
496 do l = 1, lz
497 do j = 1, lx
498 msk_c = msk_c + 1
499 this%msk(msk_c) = linear_index(j, ly, l, el, lx, ly, lz)
500 this%facet(msk_c) = 4
501 end do
502 end do
503 case (5)
504 do k = 1, ly
505 do j = 1, lx
506 msk_c = msk_c + 1
507 this%msk(msk_c) = linear_index(j, k, 1, el, lx, ly, lz)
508 this%facet(msk_c) = 5
509 end do
510 end do
511 case (6)
512 do k = 1, ly
513 do j = 1, lx
514 msk_c = msk_c + 1
515 this%msk(msk_c) = linear_index(j, k, lz, el, lx, ly, lz)
516 this%facet(msk_c) = 6
517 end do
518 end do
519 end select
520 end do
521 this%facet(0) = msk_c
522 if (neko_bcknd_device .eq. 1) then
523 !Observe the facet_mask is junk if only_facet is false
524 n = msk_c + 1
525 call device_map(this%facet, this%facet_d, n)
526 call device_memcpy(this%facet, this%facet_d, n, &
527 host_to_device, sync = .true.)
528 end if
529 if ( .not. only_facet) then
530 !Makes check for points not on facet that should have bc applied
531 call test_field%init(this%dof)
532
533 n = test_field%size()
534 test_field%x = 0.0_rp
535 !Apply this bc once
536 do i = 1, msk_c
537 test_field%x(this%msk(i),1,1,1) = 1.0
538 end do
539 if (neko_bcknd_device .eq. 1) then
540 call device_memcpy(test_field%x, test_field%x_d, n, &
541 host_to_device, sync = .true.)
542 end if
543 !Check if some point that was not zeroed was zeroed on another element
544 call this%coef%gs_h%op(test_field, gs_op_add)
545 if (neko_bcknd_device .eq. 1) then
546 call device_memcpy(test_field%x, test_field%x_d, n, &
547 device_to_host, sync = .true.)
548 end if
549 msk_c = 0
550 do i = 1, this%dof%size()
551 if (test_field%x(i,1,1,1) .gt. 0.5) then
552 msk_c = msk_c + 1
553 end if
554 end do
555 !Allocate new mask
556 deallocate(this%msk)
557 allocate(this%msk(0:msk_c))
558 j = 1
559 do i = 1, this%dof%size()
560 if (test_field%x(i,1,1,1) .gt. 0.5) then
561 this%msk(j) = i
562 j = j + 1
563 end if
564 end do
565
566 call test_field%free()
567 end if
568
569 this%msk(0) = msk_c
570
571 if (neko_bcknd_device .eq. 1) then
572 n = msk_c + 1
573 call device_map(this%msk, this%msk_d, n)
574 call device_memcpy(this%msk, this%msk_d, n, &
575 host_to_device, sync = .true.)
576 end if
577
578 if (.not. allocated(this%name)) then
579 this%name = ""
580 else
581 write(log_buf, '(A,A)') 'BC assigned name : ', trim(this%name)
582 call neko_log%message(log_buf)
583 end if
584
585 if (.not. allocated(this%zone_indices)) then
586 allocate(this%zone_indices(1))
587 this%zone_indices(1) = -1
588 end if
589
590 end subroutine bc_finalize_base
591
595 subroutine bc_debug_mask(this, file_name)
596 class(bc_t), intent(inout) :: this
597 character(len=*), intent(in) :: file_name
598 type(field_t) :: bdry_field
599 integer:: i, m, k
600 type(file_t) :: dump_file
601
602 call bdry_field%init(this%coef%dof, 'bdry')
603 m = this%msk(0)
604 do i = 1, m
605 k = this%msk(i)
606 bdry_field%x(k,1,1,1) = 1.0_rp
607 end do
608 call dump_file%init(file_name)
609 call dump_file%write(bdry_field)
610
611 end subroutine bc_debug_mask
612end 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
Defines a boundary condition.
Definition bc.f90:34
subroutine bc_mark_zone(this, bc_zone)
Mark all facets from a zone.
Definition bc.f90:418
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:359
subroutine bc_finalize_base(this, only_facets)
Finalize the construction of the bc by populting the msk and facet arrays.
Definition bc.f90:434
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:315
subroutine bc_debug_mask(this, file_name)
Write a field showing the mask of the bc.
Definition bc.f90:596
subroutine bc_mark_facet(this, facet, el)
Mark facet on element el as part of the boundary condition.
Definition bc.f90:390
subroutine bc_mark_facets(this, facet_list)
Mark all facets from a (facet, el) tuple list.
Definition bc.f90:403
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:219
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:77
integer, parameter, public log_size
Definition log.f90:46
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:35
Module with things related to the simulation time.
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: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:56
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:55
The function space for the SEM solution fields.
Definition space.f90:63
Integer 2-tuple based stack.
Definition stack.f90:84
A struct that contains all info about the time, expand as needed.
Integer based 2-tuple.
Definition tuple.f90:51