Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
bc_list.f90
Go to the documentation of this file.
1! Copyright (c) 2024-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_list
36 use num_types, only : rp
37 use field, only : field_t
39 use utils, only : neko_error
40 use, intrinsic :: iso_c_binding, only : c_ptr
41 use bc, only : bc_t, bc_ptr_t
42 use time_state, only : time_state_t
43 !$ use omp_lib
44 implicit none
45 private
46
49 type, public :: bc_list_t
50 ! The items of the list.
51 class(bc_ptr_t), allocatable, private :: items(:)
53 integer, private :: size_ = 0
56 integer, private :: capacity
57 contains
59 procedure, pass(this) :: init => bc_list_init
61 procedure, pass(this) :: free => bc_list_free
62
64 procedure, pass(this) :: append => bc_list_append
66 procedure, pass(this) :: get => bc_list_get
68 procedure, pass(this) :: get_by_name => bc_list_get_by_name
70 procedure, pass(this) :: get_by_zone_index => bc_list_get_by_zone_index
71
73 procedure, pass(this) :: is_empty => bc_list_is_empty
75 procedure, pass(this) :: strong => bc_list_strong
77 procedure :: size => bc_list_size
78
80 generic :: apply => apply_scalar, apply_vector, &
81 apply_scalar_device, apply_vector_device, &
82 apply_scalar_field, apply_vector_field
84 procedure, pass(this) :: apply_scalar => bc_list_apply_scalar_array
86 procedure, pass(this) :: apply_vector => bc_list_apply_vector_array
88 procedure, pass(this) :: apply_scalar_device => bc_list_apply_scalar_device
90 procedure, pass(this) :: apply_vector_device => bc_list_apply_vector_device
92 procedure, pass(this) :: apply_scalar_field => bc_list_apply_scalar_field
94 procedure, pass(this) :: apply_vector_field => bc_list_apply_vector_field
95 end type bc_list_t
96
97contains
98
101 subroutine bc_list_init(this, capacity)
102 class(bc_list_t), intent(inout), target :: this
103 integer, optional :: capacity
104 integer :: n
105
106 call this%free()
107
108 n = 1
109 if (present(capacity)) n = capacity
110
111 allocate(this%items(n))
112
113 this%size_ = 0
114 this%capacity = n
115
116 end subroutine bc_list_init
117
121 subroutine bc_list_free(this)
122 class(bc_list_t), intent(inout) :: this
123 integer :: i
124
125 if (allocated(this%items)) then
126 do i = 1, this%size_
127 this%items(i)%ptr => null()
128 end do
129
130 deallocate(this%items)
131 end if
132
133 this%size_ = 0
134 this%capacity = 0
135 end subroutine bc_list_free
136
140 subroutine bc_list_append(this, bc)
141 class(bc_list_t), intent(inout) :: this
142 class(bc_t), intent(inout), target :: bc
143 type(bc_ptr_t), allocatable :: tmp(:)
144
145 if (this%size_ .ge. this%capacity) then
146 this%capacity = this%capacity * 2
147 allocate(tmp(this%capacity))
148 tmp(1:this%size_) = this%items
149 call move_alloc(tmp, this%items)
150 end if
151
152 this%size_ = this%size_ + 1
153 this%items(this%size_)%ptr => bc
154
155 end subroutine bc_list_append
156
160 function bc_list_get(this, i) result(bc)
161 class(bc_list_t), intent(in) :: this
162 class(bc_t), pointer :: bc
163 integer, intent(in) :: i
164
165 if (i .lt. 1 .or. i .gt. this%size_) then
166 call neko_error("Index out of bounds in bc_list_get")
167 end if
168
169 bc => this%items(i)%ptr
170
171 end function bc_list_get
172
176 function bc_list_get_by_name(this, name) result(bc)
177 class(bc_list_t), intent(in) :: this
178 class(bc_t), pointer :: bc
179 character(len=*), intent(in) :: name
180 integer :: i
181
182 do i = 1, this%size_
183 if (this%items(i)%ptr%name .eq. trim(name)) then
184 bc => this%items(i)%ptr
185 return
186 end if
187 end do
188
189 ! If the function reaches this point, no item was found
190 call neko_error("Name not found in bc_list")
191
192 end function bc_list_get_by_name
193
197 function bc_list_get_by_zone_index(this, zone_index) result(bc)
198 class(bc_list_t), intent(in) :: this
199 class(bc_t), pointer :: bc
200 integer, intent(in) :: zone_index
201 integer :: i, j
202
203 do i = 1, this%size_
204 do j = 1, size(this%items(i)%ptr%zone_indices)
205 if (this%items(i)%ptr%zone_indices(j) == zone_index) then
206 bc => this%items(i)%ptr
207 return
208 end if
209 end do
210 end do
211
212 ! If the function reaches this point, no item was found
213 call neko_error("Zone index not found in bc_list")
214
215 end function bc_list_get_by_zone_index
216
224 subroutine bc_list_apply_scalar_array(this, x, n, time, strong, strm)
225 class(bc_list_t), intent(inout) :: this
226 integer, intent(in) :: n
227 real(kind=rp), intent(inout), dimension(n) :: x
228 type(time_state_t), intent(in), optional :: time
229 logical, intent(in), optional :: strong
230 type(c_ptr), intent(inout), optional :: strm
231 logical :: strong_
232 type(c_ptr) :: x_d
233 integer :: i
234
235 if (neko_bcknd_device .eq. 1) then
236
237 x_d = device_get_ptr(x)
238
239 call this%apply_scalar_device(x_d, time = time, &
240 strong = strong, strm = strm)
241 else
242 ! Resolve strong into a concrete, always-present local before opening
243 ! the parallel region. CCE's outlined region prologue dereferences a
244 ! null descriptor for any absent optional dummy referenced inside the
245 ! region, and it does not treat an unallocated allocatable as not
246 ! present, so optionals must not be forwarded into the region. time
247 ! cannot be given a meaningful concrete default, so we branch on
248 ! present(time) outside the region instead.
249 strong_ = .true.
250 if (present(strong)) strong_ = strong
251
252 if (present(time)) then
253 !$omp parallel
254 do i = 1, this%size_
255 call this%items(i)%ptr%apply_scalar(x, n, time = time, &
256 strong = strong_)
257 end do
258 !$omp end parallel
259 else
260 !$omp parallel
261 do i = 1, this%size_
262 call this%items(i)%ptr%apply_scalar(x, n, strong = strong_)
263 end do
264 !$omp end parallel
265 end if
266 end if
267 end subroutine bc_list_apply_scalar_array
268
278 subroutine bc_list_apply_vector_array(this, x, y, z, n, time, strong, strm)
279 class(bc_list_t), intent(inout) :: this
280 integer, intent(in) :: n
281 real(kind=rp), intent(inout), dimension(n) :: x
282 real(kind=rp), intent(inout), dimension(n) :: y
283 real(kind=rp), intent(inout), dimension(n) :: z
284 type(time_state_t), intent(in), optional :: time
285 logical, intent(in), optional :: strong
286 type(c_ptr), intent(inout), optional :: strm
287 logical :: strong_
288 type(c_ptr) :: x_d
289 type(c_ptr) :: y_d
290 type(c_ptr) :: z_d
291 integer :: i
292
293 if (neko_bcknd_device .eq. 1) then
294
295 x_d = device_get_ptr(x)
296 y_d = device_get_ptr(y)
297 z_d = device_get_ptr(z)
298
299 call this%apply_vector_device(x_d, y_d, z_d, time = time, &
300 strong = strong, strm = strm)
301 else
302 ! Resolve strong into a concrete, always-present local before opening
303 ! the parallel region. CCE's outlined region prologue dereferences a
304 ! null descriptor for any absent optional dummy referenced inside the
305 ! region, and it does not treat an unallocated allocatable as not
306 ! present, so optionals must not be forwarded into the region. time
307 ! cannot be given a meaningful concrete default, so we branch on
308 ! present(time) outside the region instead.
309 strong_ = .true.
310 if (present(strong)) strong_ = strong
311
312 if (present(time)) then
313 !$omp parallel
314 do i = 1, this%size_
315 call this%items(i)%ptr%apply_vector(x, y, z, n, time = time, &
316 strong = strong_)
317 end do
318 !$omp end parallel
319 else
320 !$omp parallel
321 do i = 1, this%size_
322 call this%items(i)%ptr%apply_vector(x, y, z, n, strong = strong_)
323 end do
324 !$omp end parallel
325 end if
326 end if
327
328 end subroutine bc_list_apply_vector_array
329
336 subroutine bc_list_apply_scalar_device(this, x_d, time, strong, strm)
337 class(bc_list_t), intent(inout) :: this
338 type(c_ptr), intent(inout) :: x_d
339 type(time_state_t), intent(in), optional :: time
340 logical, intent(in), optional :: strong
341 type(c_ptr), intent(inout), optional :: strm
342 type(c_ptr) :: strm_
343 integer :: i
344
345 if (present(strm)) then
346 strm_ = strm
347 else
348 strm_ = glb_cmd_queue
349 end if
350
351 do i = 1, this%size_
352 call this%items(i)%ptr%apply_scalar_dev(x_d, time = time, &
353 strong = strong, strm = strm_)
354 end do
355
356 end subroutine bc_list_apply_scalar_device
357
366 subroutine bc_list_apply_vector_device(this, x_d, y_d, z_d, time, strong, &
367 strm)
368 class(bc_list_t), intent(inout) :: this
369 type(c_ptr), intent(inout) :: x_d
370 type(c_ptr), intent(inout) :: y_d
371 type(c_ptr), intent(inout) :: z_d
372 type(time_state_t), intent(in), optional :: time
373 logical, intent(in), optional :: strong
374 type(c_ptr), intent(inout), optional :: strm
375 type(c_ptr) :: strm_
376 integer :: i
377
378 if (present(strm)) then
379 strm_ = strm
380 else
381 strm_ = glb_cmd_queue
382 end if
383
384 do i = 1, this%size_
385 call this%items(i)%ptr%apply_vector_dev(x_d, y_d, z_d, time = time, &
386 strong = strong, strm = strm_)
387 end do
388
389 end subroutine bc_list_apply_vector_device
390
397 subroutine bc_list_apply_scalar_field(this, x, time, strong, strm)
398 class(bc_list_t), intent(inout) :: this
399 type(field_t), intent(inout) :: x
400 type(time_state_t), intent(in), optional :: time
401 logical, intent(in), optional :: strong
402 type(c_ptr), intent(inout), optional :: strm
403 logical :: strong_
404 type(c_ptr) :: strm_
405 integer :: i
406
407 ! Resolve strong and strm into concrete, always-present locals before
408 ! opening the parallel region. CCE's outlined region prologue dereferences
409 ! a null descriptor for any absent optional dummy referenced inside the
410 ! region, and it does not treat an unallocated allocatable as not present,
411 ! so optionals must not be forwarded into the region. time cannot be given
412 ! a meaningful concrete default, so we branch on present(time) outside the
413 ! region instead.
414 strong_ = .true.
415 if (present(strong)) strong_ = strong
416 strm_ = glb_cmd_queue
417 if (present(strm)) strm_ = strm
418
419 if (present(time)) then
420 !$omp parallel if (.not. omp_in_parallel())
421 do i = 1, this%size_
422 call this%items(i)%ptr%apply_scalar_generic(x, time = time, &
423 strong = strong_, strm = strm_)
424 end do
425 !$omp end parallel
426 else
427 !$omp parallel if (.not. omp_in_parallel())
428 do i = 1, this%size_
429 call this%items(i)%ptr%apply_scalar_generic(x, &
430 strong = strong_, strm = strm_)
431 end do
432 !$omp end parallel
433 end if
434
435 end subroutine bc_list_apply_scalar_field
436
445 subroutine bc_list_apply_vector_field(this, x, y, z, time, strong, strm)
446 class(bc_list_t), intent(inout) :: this
447 type(field_t), intent(inout) :: x
448 type(field_t), intent(inout) :: y
449 type(field_t), intent(inout) :: z
450 type(time_state_t), intent(in), optional :: time
451 logical, intent(in), optional :: strong
452 type(c_ptr), intent(inout), optional :: strm
453 logical :: strong_
454 type(c_ptr) :: strm_
455 integer :: i
456
457 ! Resolve strong and strm into concrete, always-present locals before
458 ! opening the parallel region. CCE's outlined region prologue dereferences
459 ! a null descriptor for any absent optional dummy referenced inside the
460 ! region, and it does not treat an unallocated allocatable as not present,
461 ! so optionals must not be forwarded into the region. time cannot be given
462 ! a meaningful concrete default, so we branch on present(time) outside the
463 ! region instead.
464 strong_ = .true.
465 if (present(strong)) strong_ = strong
466 strm_ = glb_cmd_queue
467 if (present(strm)) strm_ = strm
468
469 if (present(time)) then
470 !$omp parallel if (.not. omp_in_parallel())
471 do i = 1, this%size_
472 call this%items(i)%ptr%apply_vector_generic(x, y, z, time = time, &
473 strong = strong_, strm = strm_)
474 end do
475 !$omp end parallel
476 else
477 !$omp parallel if (.not. omp_in_parallel())
478 do i = 1, this%size_
479 call this%items(i)%ptr%apply_vector_generic(x, y, z, &
480 strong = strong_, strm = strm_)
481 end do
482 !$omp end parallel
483 end if
484
485 end subroutine bc_list_apply_vector_field
486
488 pure function bc_list_strong(this, i) result(strong)
489 class(bc_list_t), intent(in), target :: this
490 integer, intent(in) :: i
491 logical :: strong
492
493 strong = this%items(i)%ptr%strong
494 end function bc_list_strong
495
497 function bc_list_is_empty(this) result(is_empty)
498 class(bc_list_t), intent(in), target :: this
499 logical :: is_empty
500 integer :: i
501
502 is_empty = .true.
503 do i = 1, this%size_
504
505 if (.not. allocated(this%items(i)%ptr%msk)) then
506 call neko_error("bc not finalized, error in bc_list%is_empty")
507 end if
508
509 if (this%items(i)%ptr%msk(0) > 0) is_empty = .false.
510
511 end do
512 end function bc_list_is_empty
513
515 pure function bc_list_size(this) result(size)
516 class(bc_list_t), intent(in), target :: this
517 integer :: size
518
519 size = this%size_
520 end function bc_list_size
521
522end module bc_list
Return the device pointer for an associated Fortran array.
Definition device.F90:108
Defines a list of bc_t.
Definition bc_list.f90:34
pure logical function bc_list_strong(this, i)
Return whether the bc is strong or not.
Definition bc_list.f90:489
subroutine bc_list_apply_vector_array(this, x, y, z, n, time, strong, strm)
Apply a list of boundary conditions to a vector field.
Definition bc_list.f90:279
subroutine bc_list_apply_vector_device(this, x_d, y_d, z_d, time, strong, strm)
Apply a list of boundary conditions to a vector field on the device.
Definition bc_list.f90:368
class(bc_t) function, pointer bc_list_get(this, i)
Get the item at the given index.
Definition bc_list.f90:161
subroutine bc_list_apply_scalar_array(this, x, n, time, strong, strm)
Apply a list of boundary conditions to a scalar field.
Definition bc_list.f90:225
subroutine bc_list_apply_scalar_device(this, x_d, time, strong, strm)
Apply a list of boundary conditions to a scalar field on the device.
Definition bc_list.f90:337
pure integer function bc_list_size(this)
Return the number of items in the list.
Definition bc_list.f90:516
subroutine bc_list_append(this, bc)
Append a condition to the end of the list.
Definition bc_list.f90:141
subroutine bc_list_init(this, capacity)
Constructor.
Definition bc_list.f90:102
class(bc_t) function, pointer bc_list_get_by_name(this, name)
Get the item from a given name.
Definition bc_list.f90:177
subroutine bc_list_apply_vector_field(this, x, y, z, time, strong, strm)
Apply a list of boundary conditions to a vector field.
Definition bc_list.f90:446
subroutine bc_list_free(this)
Destructor.
Definition bc_list.f90:122
class(bc_t) function, pointer bc_list_get_by_zone_index(this, zone_index)
Get the item from zone_index.
Definition bc_list.f90:198
subroutine bc_list_apply_scalar_field(this, x, time, strong, strm)
Apply a list of boundary conditions to a scalar field.
Definition bc_list.f90:398
logical function bc_list_is_empty(this)
Return whether the list is empty.
Definition bc_list.f90:498
Defines a boundary condition.
Definition bc.f90:34
Device abstraction, common interface for various accelerators.
Definition device.F90:34
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:52
Defines a field.
Definition field.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
Pointer to a `bc_t`.
Definition bc.f90:133
Base type for a boundary condition.
Definition bc.f90:62
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:49
A struct that contains all info about the time, expand as needed.