Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
overset_interface_vector.f90
Go to the documentation of this file.
1! Copyright (c) 2026, 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!
35 use comm, only : neko_global_comm
37 use registry, only : neko_registry
38 use num_types, only : rp
39 use coefs, only : coef_t
40 use dirichlet, only : dirichlet_t
43 use mask, only : mask_t
44 use dofmap, only : dofmap_t
45 use bc, only : bc_t
46 use bc_list, only : bc_list_t
47 use utils, only : split_string
48 use field, only : field_t
49 use field_list, only : field_list_t
50 use math, only : masked_copy_0
52 use dofmap, only : dofmap_t
53 use vector, only : vector_t
54 use vector_list, only : vector_list_t
56 use math, only : copy
58 use vector_math, only : vector_copy
62 use stack, only: stack_i4_t
63 use json_module, only : json_file
65 use field_list, only : field_list_t
66 use, intrinsic :: iso_c_binding, only : c_ptr, c_size_t
67 use time_state, only : time_state_t
68 implicit none
69 private
70
72 ! for the application on a vector field.
73 type, public, extends(bc_t) :: overset_interface_vector_t
74 ! The bc for the first compoent.
75 type(field_dirichlet_t) :: bc_u
76 ! The bc for the second compoent.
77 type(field_dirichlet_t) :: bc_v
78 ! The bc for the third compoent.
79 type(field_dirichlet_t) :: bc_w
83 type(global_interpolation_t) :: interface_interpolator
85 type(mask_t) :: interface_dof_mask
86 type(mask_t) :: domain_element_mask
88 type(vector_t) :: x_dof, y_dof, z_dof
89 type(vector_t) :: x_interface_dof, y_interface_dof, z_interface_dof
90 type(vector_t) :: u_interface, v_interface, w_interface
91 type(vector_list_t) :: interface_dof, interface_field
93 type(global_interpolation_settings_t) :: interpolation_settings
94 logical :: find_interface = .false.
95 logical :: setup = .false.
96
99 procedure(morph_overset_interface), nopass, pointer :: morph_interface => null()
100
101 contains
103 procedure, pass(this) :: init => overset_interface_vector_init
105 procedure, pass(this) :: init_from_components => &
108 procedure, pass(this) :: free => overset_interface_vector_free
110 procedure, pass(this) :: finalize => overset_interface_vector_finalize
112 procedure, pass(this) :: apply_scalar => &
115 procedure, pass(this) :: apply_vector => &
118 procedure, pass(this) :: apply_vector_dev => &
121 procedure, pass(this) :: apply_scalar_dev => &
123 procedure, pass(this) :: update => overset_interface_update
124
126 procedure, pass(this), private :: build_masks_ => build_masks_
128 procedure, pass(this), private :: gather_interface_dofs_ => gather_interface_dofs_
130 procedure, pass(this), private :: setup_interpolator_ => setup_interpolator_
131
133
134contains
135
139 subroutine overset_interface_vector_init(this, coef, json)
140 class(overset_interface_vector_t), intent(inout), target :: this
141 type(coef_t), target, intent(in) :: coef
142 type(json_file), intent(inout) ::json
143 real(kind=rp) :: tol, pad
144
146 call json_get_or_default(json, "interpolation.tolerance", &
147 tol, -1.0_rp)
148 call json_get_or_default(json, "interpolation.padding", &
149 pad, -1.0_rp)
150
151 call this%init_from_components(coef, tol, pad)
152
153 end subroutine overset_interface_vector_init
154
157 subroutine overset_interface_vector_init_from_components(this, coef, tol, pad)
158 class(overset_interface_vector_t), intent(inout), target :: this
159 type(coef_t), intent(in) :: coef
160 real(kind=rp), intent(in), optional :: tol, pad
161
163 call this%init_base(coef)
164
166 if (present(tol)) then
167 if (tol .gt. 0.0_rp) this%interpolation_settings%tolerance = tol
168 end if
169 if (present(pad)) then
170 if (pad .gt. 0.0_rp) this%interpolation_settings%padding = pad
171 end if
172
173 call this%bc_u%init_from_components(coef, "u")
174 call this%bc_v%init_from_components(coef, "v")
175 call this%bc_w%init_from_components(coef, "w")
176
177 call this%field_list%init(3)
178 call this%field_list%assign_to_field(1, this%bc_u%field_bc)
179 call this%field_list%assign_to_field(2, this%bc_v%field_bc)
180 call this%field_list%assign_to_field(3, this%bc_w%field_bc)
181
183 call this%x_dof%init(this%dof%size(), 'x')
184 call this%y_dof%init(this%dof%size(), 'y')
185 call this%z_dof%init(this%dof%size(), 'z')
186
189 if (neko_bcknd_device .eq. 1) then
190 ! copy
191 call device_copy(this%x_dof%x_d, this%dof%x_d, this%dof%size())
192 call device_copy(this%y_dof%x_d, this%dof%y_d, this%dof%size())
193 call device_copy(this%z_dof%x_d, this%dof%z_d, this%dof%size())
194 ! synchronize
195 call this%x_dof%copy_from(device_to_host, sync = .false.)
196 call this%y_dof%copy_from(device_to_host, sync = .false.)
197 call this%z_dof%copy_from(device_to_host, sync = .true.)
198 else
199 call copy(this%x_dof%x, this%dof%x, this%dof%size())
200 call copy(this%y_dof%x, this%dof%y, this%dof%size())
201 call copy(this%z_dof%x, this%dof%z, this%dof%size())
202 end if
203
205
209 class(overset_interface_vector_t), target, intent(inout) :: this
210
211 call this%bc_u%free()
212 call this%bc_v%free()
213 call this%bc_w%free()
214
215 call this%field_list%free()
216 call this%interface_dof%free()
217 call this%interface_field%free()
218
219 call this%x_dof%free()
220 call this%y_dof%free()
221 call this%z_dof%free()
222
223 call this%x_interface_dof%free()
224 call this%y_interface_dof%free()
225 call this%z_interface_dof%free()
226 call this%u_interface%free()
227 call this%v_interface%free()
228 call this%w_interface%free()
229 call this%interface_interpolator%free()
230 call this%interface_dof_mask%free()
231 call this%domain_element_mask%free()
232 call this%free_base()
233
234 !if (associated(this%update_)) then
235 ! nullify(this%update_)
236 !end if
237 end subroutine overset_interface_vector_free
238
243 subroutine overset_interface_vector_apply_scalar(this, x, n, time, strong)
244 class(overset_interface_vector_t), intent(inout) :: this
245 integer, intent(in) :: n
246 real(kind=rp), intent(inout), dimension(n) :: x
247 type(time_state_t), intent(in), optional :: time
248 logical, intent(in), optional :: strong
249
250 call neko_error("overset_interface_vector cannot apply scalar BCs.&
251 & Use overset_interface_vector::apply_vector instead!")
252
254
258 subroutine overset_interface_vector_apply_scalar_dev(this, x_d, time, &
259 strong, strm)
260 class(overset_interface_vector_t), intent(inout), target :: this
261 type(c_ptr), intent(inout) :: x_d
262 type(time_state_t), intent(in), optional :: time
263 logical, intent(in), optional :: strong
264 type(c_ptr), intent(inout) :: strm
265
266 call neko_error("overset_interface_vector cannot apply scalar BCs.&
267 & Use overset_interface_vector::apply_vector instead!")
268
270
277 subroutine overset_interface_vector_apply_vector(this, x, y, z, n, time, &
278 strong)
279 class(overset_interface_vector_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 logical :: strong_
287
288 if (present(strong)) then
289 strong_ = strong
290 else
291 strong_ = .true.
292 end if
293
294 if (strong_) then
295
296 ! We can send any of the 3 bcs we have as argument, since they are all
297 ! the same boundary.
298 if (.not. this%updated) then
299 call this%update(time)
300 this%updated = .true.
301 end if
302
304 call masked_copy_0(x, this%bc_u%field_bc%x, this%msk, n, this%msk(0))
305 call masked_copy_0(y, this%bc_v%field_bc%x, this%msk, n, this%msk(0))
306 call masked_copy_0(z, this%bc_w%field_bc%x, this%msk, n, this%msk(0))
307 end if
308
310
317 subroutine overset_interface_vector_apply_vector_dev(this, x_d, y_d, z_d, &
318 time, strong, strm)
319 class(overset_interface_vector_t), intent(inout), target :: this
320 type(c_ptr), intent(inout) :: x_d
321 type(c_ptr), intent(inout) :: y_d
322 type(c_ptr), intent(inout) :: z_d
323 type(time_state_t), intent(in), optional :: time
324 logical, intent(in), optional :: strong
325 type(c_ptr), intent(inout) :: strm
326 logical :: strong_
327
328 if (present(strong)) then
329 strong_ = strong
330 else
331 strong_ = .true.
332 end if
333
334 if (strong_) then
335 if (.not. this%updated) then
336 call this%update(time)
337 this%updated = .true.
338 end if
339
340 if (this%msk(0) .gt. 0) then
341 call device_masked_copy_0(x_d, this%bc_u%field_bc%x_d, &
342 this%bc_u%msk_d, this%bc_u%dof%size(), this%msk(0), strm) ! adperez: change the masks used here
343 call device_masked_copy_0(y_d, this%bc_v%field_bc%x_d, &
344 this%bc_v%msk_d, this%bc_v%dof%size(), this%msk(0), strm)
345 call device_masked_copy_0(z_d, this%bc_w%field_bc%x_d, &
346 this%bc_w%msk_d, this%bc_w%dof%size(), this%msk(0), strm)
347 end if
348 end if
349
351
353 subroutine overset_interface_vector_finalize(this, only_facets)
354 class(overset_interface_vector_t), target, intent(inout) :: this
355 logical, optional, intent(in) :: only_facets
356 logical :: only_facets_
357
358 if (present(only_facets)) then
359 only_facets_ = only_facets
360 else
361 only_facets_ = .false.
362 end if
363
365 call this%finalize_base(only_facets_)
366
367 call this%bc_u%mark_facets(this%marked_facet)
368 call this%bc_v%mark_facets(this%marked_facet)
369 call this%bc_w%mark_facets(this%marked_facet)
370
371 call this%bc_u%finalize(only_facets_)
372 call this%bc_v%finalize(only_facets_)
373 call this%bc_w%finalize(only_facets_)
374
376 call this%build_masks_()
377
379 call this%x_interface_dof%init(this%interface_dof_mask%size(), 'x_interface')
380 call this%y_interface_dof%init(this%interface_dof_mask%size(), 'y_interface')
381 call this%z_interface_dof%init(this%interface_dof_mask%size(), 'z_interface')
382 call this%gather_interface_dofs_()
383
385 call this%setup_interpolator_()
386
388 call this%u_interface%init(this%interface_dof_mask%size(), 'u_interface')
389 call this%v_interface%init(this%interface_dof_mask%size(), 'v_interface')
390 call this%w_interface%init(this%interface_dof_mask%size(), 'w_interface')
391
393 call this%interface_dof%init(3)
394 call this%interface_dof%assign_to_vector(1, this%x_interface_dof)
395 call this%interface_dof%assign_to_vector(2, this%y_interface_dof)
396 call this%interface_dof%assign_to_vector(3, this%z_interface_dof)
397
398 call this%interface_field%init(3)
399 call this%interface_field%assign_to_vector(1, this%u_interface)
400 call this%interface_field%assign_to_vector(2, this%v_interface)
401 call this%interface_field%assign_to_vector(3, this%w_interface)
402
403
405
407 subroutine overset_interface_update(this, time)
408 class(overset_interface_vector_t), intent(inout) :: this
409 type(time_state_t), intent(in) :: time
410 type(field_t), pointer :: u, v, w
411
412
414 call this%morph_interface(this%interface_dof, this%interface_field, &
415 this%interface_dof_mask, time, this%name, &
416 this%find_interface)
417
419 ! not implemented for now
420 ! if (substep .eq. 1) then
421 ! call this%extrapolate()
422 ! end if
423
425 if (this%find_interface) then
426
427 ! sync
428 call this%x_interface_dof%copy_from(device_to_host, sync = .false.)
429 call this%y_interface_dof%copy_from(device_to_host, sync = .false.)
430 call this%z_interface_dof%copy_from(device_to_host, sync = .true.)
431
432 call this%interface_interpolator%find_points(this%x_interface_dof%x, &
433 this%y_interface_dof%x, this%z_interface_dof%x, &
434 this%x_interface_dof%size())
435 this%find_interface = .false.
436
437 end if
438
440 u => neko_registry%get_field("u")
441 v => neko_registry%get_field("v")
442 w => neko_registry%get_field("w")
443
445 call this%interface_interpolator%evaluate_masked(this%u_interface%x, u%x, this%domain_element_mask, .false.)
446 call this%interface_interpolator%evaluate_masked(this%v_interface%x, v%x, this%domain_element_mask, .false.)
447 call this%interface_interpolator%evaluate_masked(this%w_interface%x, w%x, this%domain_element_mask, .false.)
448
450 call vector_masked_scatter_copy(this%bc_u%field_bc%x(:,1,1,1), this%u_interface, &
451 this%interface_dof_mask, this%bc_u%dof%size())
452 call vector_masked_scatter_copy(this%bc_v%field_bc%x(:,1,1,1), this%v_interface, &
453 this%interface_dof_mask, this%bc_v%dof%size())
454 call vector_masked_scatter_copy(this%bc_w%field_bc%x(:,1,1,1), this%w_interface, &
455 this%interface_dof_mask, this%bc_w%dof%size())
456
457
458 end subroutine overset_interface_update
459
460 !===================
461 ! Helper subroutines
462 !===================
463
465 subroutine build_masks_(this)
466 class(overset_interface_vector_t), intent(inout) :: this
467 type(mask_t) :: temp_mask
468 logical, allocatable :: found(:)
469 integer :: i, j, k, e, new_size, nelems
470 integer :: lx, ly, lz
471 integer :: nonlinear_idx(4), linear_idx
472 type(stack_i4_t) :: stack
473
475 call this%interface_dof_mask%init(this%msk(1:this%msk(0)), this%msk(0))
476
479 lx = this%Xh%lx
480 ly = this%Xh%ly
481 lz = this%Xh%lz
482 allocate(found(this%msh%nelv))
483 found = .false.
484 !! Find sem elements that contain the boundary points
485 do i = 1, this%msk(0)
486 linear_idx = this%msk(i)
487 nonlinear_idx = nonlinear_index(linear_idx, lx, ly, lz)
488 found(nonlinear_idx(4)) = .true.
489 end do
490 !! fill the stack containing the gll indices
491 nelems = 0
492 call stack%init()
493 do e = 1, this%msh%nelv
494 if (found(e)) then
495 nelems = nelems + 1
496 do k = 1, this%Xh%lz
497 do j = 1, this%Xh%ly
498 do i = 1, this%Xh%lx
499 linear_idx = linear_index(i, j, k, e, lx, ly, lz)
500 call stack%push(linear_idx)
501 end do
502 end do
503 end do
504 end if
505 end do
506 deallocate(found)
507 call temp_mask%init(stack%array(), stack%size())
508 call stack%free()
509
511 call this%domain_element_mask%invert_mask(temp_mask, this%dof%size())
512
513 call temp_mask%free()
514
515 end subroutine build_masks_
516
518 subroutine gather_interface_dofs_(this)
519 class(overset_interface_vector_t), intent(inout) :: this
520
522 call vector_masked_gather_copy(this%x_interface_dof, &
523 this%dof%x(:,1,1,1), &
524 this%interface_dof_mask, &
525 this%dof%size())
526 call vector_masked_gather_copy(this%y_interface_dof, &
527 this%dof%y(:,1,1,1), &
528 this%interface_dof_mask, &
529 this%dof%size())
530 call vector_masked_gather_copy(this%z_interface_dof, &
531 this%dof%z(:,1,1,1), &
532 this%interface_dof_mask, &
533 this%dof%size())
535 call this%x_interface_dof%copy_from(device_to_host, sync = .false.)
536 call this%y_interface_dof%copy_from(device_to_host, sync = .false.)
537 call this%z_interface_dof%copy_from(device_to_host, sync = .true.)
538
539 end subroutine gather_interface_dofs_
540
541
543 subroutine setup_interpolator_(this)
544 class(overset_interface_vector_t), intent(inout) :: this
545
547 call this%interface_interpolator%init(this%dof, &
549 tol=this%interpolation_settings%tolerance, &
550 pad=this%interpolation_settings%padding, &
551 mask=this%domain_element_mask)
552
554 call this%interface_interpolator%find_points(this%x_interface_dof%x, &
555 this%y_interface_dof%x, &
556 this%z_interface_dof%x, &
557 this%x_interface_dof%size())
558
559 end subroutine setup_interpolator_
560
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
Abstract interface defining a dirichlet condition on a list of fields.
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
User callback for overset-interface morphing and boundary-value updates.
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_global_comm
Definition comm.F90:44
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_masked_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Copy a masked vector .
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
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines user dirichlet condition for a scalar field.
Defines a field.
Definition field.f90:34
Implements global_interpolation given a dofmap.
Utilities for retrieving parameters from the case files.
Object for handling masks in Neko.
Definition mask.f90:34
Definition math.f90:60
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:252
subroutine, public masked_copy_0(a, b, mask, n, n_mask)
Copy a masked vector .
Definition math.f90:271
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines overset interface vector boundary conditions.
subroutine overset_interface_vector_finalize(this, only_facets)
Finalize by building the mask arrays and propagating to underlying bcs.
subroutine overset_interface_vector_init_from_components(this, coef, tol, pad)
Constructor from components.
subroutine overset_interface_vector_apply_scalar_dev(this, x_d, time, strong, strm)
No-op apply scalar (device).
subroutine overset_interface_vector_free(this)
Destructor. Currently unused as is, all field_dirichlet attributes are freed in fluid_scheme_incompre...
subroutine overset_interface_vector_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
Apply the boundary condition to a vector field on the device.
subroutine overset_interface_vector_apply_scalar(this, x, n, time, strong)
No-op apply scalar.
subroutine overset_interface_vector_init(this, coef, json)
Constructor.
subroutine overset_interface_vector_apply_vector(this, x, y, z, n, time, strong)
Apply the boundary condition to a vector field.
Defines overset interface scalar boundary conditions.
subroutine overset_interface_update(this, time)
Update values at the overset interface.
subroutine gather_interface_dofs_(this)
Gather interface dofs.
subroutine setup_interpolator_(this)
Set up the global interpolator.
subroutine build_masks_(this)
Build masks.
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Implements a dynamic stack ADT.
Definition stack.f90:49
Module with things related to the simulation time.
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:245
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:283
subroutine, public vector_copy(a, b, n)
Copy a vector .
subroutine, public vector_masked_gather_copy(a, b, mask, n)
Gather a vector to reduced contigous array .
subroutine, public vector_masked_scatter_copy(a, b, mask, n)
Scatter a contiguous vector into an array .
Defines a vector.
Definition vector.f90:34
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:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
User defined dirichlet condition, for which the user can work with an entire field....
field_list_t, To be able to group fields together
Implements the settings helper data container for global interpolation.
Implements global interpolation for arbitrary points in the domain.
Type for consistently handling masks in Neko. This type encapsulates the mask array and its associate...
Definition mask.f90:51
Extension of the user defined dirichlet condition overset_interface
Integer based stack.
Definition stack.f90:77
A struct that contains all info about the time, expand as needed.
vector_list_t, To be able to group vectors together