Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
overset_interface.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
42 use mask, only : mask_t
43 use bc, only : bc_t
44 use field_list, only : field_list_t
45 use math, only : masked_copy_0, copy
47 use vector, only : vector_t
48 use vector_list, only : vector_list_t
50 use device, only : device_to_host
53 use stack, only : stack_i4_t
54 use json_module, only : json_file
56 use field, only : field_t
57 use logger, only : neko_log, log_size
58 use, intrinsic :: iso_c_binding, only : c_ptr
59 use time_state, only : time_state_t
60 implicit none
61 private
62
64 type, public, extends(bc_t) :: overset_interface_t
66 type(field_dirichlet_t) :: bc_s
70 character(len=:), allocatable :: field_name
72 type(global_interpolation_t) :: interface_interpolator
74 type(mask_t) :: interface_dof_mask
75 type(mask_t) :: domain_element_mask
77 type(vector_t) :: x_dof, y_dof, z_dof
78 type(vector_t) :: x_interface_dof, y_interface_dof, z_interface_dof
80 type(vector_t) :: s_interface
81 type(vector_list_t) :: interface_dof, interface_field
83 type(global_interpolation_settings_t) :: interpolation_settings
84 logical :: find_interface = .false.
85 logical :: setup = .false.
86
89 procedure(morph_overset_interface), nopass, pointer :: morph_interface => null()
90
91 contains
93 procedure, pass(this) :: init => overset_interface_init
95 procedure, pass(this) :: init_from_components => &
98 procedure, pass(this) :: free => overset_interface_free
100 procedure, pass(this) :: finalize => overset_interface_finalize
102 procedure, pass(this) :: apply_scalar => overset_interface_apply_scalar
104 procedure, pass(this) :: apply_vector => overset_interface_apply_vector
106 procedure, pass(this) :: apply_vector_dev => overset_interface_apply_vector_dev
108 procedure, pass(this) :: apply_scalar_dev => overset_interface_apply_scalar_dev
109 procedure, pass(this) :: update => overset_interface_update
110
112 procedure, pass(this), private :: build_masks_ => build_masks_
114 procedure, pass(this), private :: gather_interface_dofs_ => gather_interface_dofs_
116 procedure, pass(this), private :: setup_interpolator_ => setup_interpolator_
117 end type overset_interface_t
118
119 abstract interface
120
135 subroutine morph_overset_interface(interface_dof, interface_field, &
136 interface_mask, time, bc_name, &
137 find_interface)
139 type(vector_list_t), intent(inout) :: interface_dof
140 type(vector_list_t), intent(inout) :: interface_field
141 type(mask_t), intent(in) :: interface_mask
142 type(time_state_t), intent(in) :: time
143 character(len=*), intent(in) :: bc_name
144 logical, intent(inout) :: find_interface
145 end subroutine morph_overset_interface
146 end interface
147
149
150contains
151
155 subroutine overset_interface_init(this, coef, json)
156 class(overset_interface_t), intent(inout), target :: this
157 type(coef_t), target, intent(in) :: coef
158 type(json_file), intent(inout) :: json
159 character(len=:), allocatable :: field_name
160 real(kind=rp) :: tol, pad
161
162 call json_get(json, "field_name", field_name)
163 call json_get_or_default(json, "interpolation.tolerance", tol, -1.0_rp)
164 call json_get_or_default(json, "interpolation.padding", pad, -1.0_rp)
165 call this%init_from_components(coef, field_name, tol, pad)
166 if (allocated(field_name)) deallocate(field_name)
167
168 end subroutine overset_interface_init
169
172 subroutine overset_interface_init_from_components(this, coef, field_name, tol, pad)
173 class(overset_interface_t), intent(inout), target :: this
174 type(coef_t), intent(in) :: coef
175 character(len=*), intent(in) :: field_name
176 real(kind=rp), intent(in), optional :: tol, pad
177 character(len=LOG_SIZE) :: log_buf
178
179 call this%init_base(coef)
180
181 if (present(tol)) then
182 if (tol .gt. 0.0_rp) then
183 this%interpolation_settings%tolerance = tol
184 end if
185 end if
186
187 if (present(pad)) then
188 if (pad .gt. 0.0_rp) then
189 this%interpolation_settings%padding = pad
190 end if
191 end if
192
193 this%field_name = field_name
194 write (log_buf, '(A,A)') "Coupling overset interface for: ", trim(this%field_name)
195 call neko_log%message(log_buf)
196
197 call this%bc_s%init_from_components(coef, this%field_name)
198 call this%field_list%init(1)
199 call this%field_list%assign_to_field(1, this%bc_s%field_bc)
200
201 call this%x_dof%init(this%dof%size(), 'x')
202 call this%y_dof%init(this%dof%size(), 'y')
203 call this%z_dof%init(this%dof%size(), 'z')
204
205 if (neko_bcknd_device .eq. 1) then
206 call device_copy(this%x_dof%x_d, this%dof%x_d, this%dof%size())
207 call device_copy(this%y_dof%x_d, this%dof%y_d, this%dof%size())
208 call device_copy(this%z_dof%x_d, this%dof%z_d, this%dof%size())
209
210 call this%x_dof%copy_from(device_to_host, sync = .false.)
211 call this%y_dof%copy_from(device_to_host, sync = .false.)
212 call this%z_dof%copy_from(device_to_host, sync = .true.)
213 else
214 call copy(this%x_dof%x, this%dof%x, this%dof%size())
215 call copy(this%y_dof%x, this%dof%y, this%dof%size())
216 call copy(this%z_dof%x, this%dof%z, this%dof%size())
217 end if
218
220
222 subroutine overset_interface_free(this)
223 class(overset_interface_t), target, intent(inout) :: this
224
225 call this%bc_s%free()
226 call this%field_list%free()
227 call this%interface_dof%free()
228 call this%interface_field%free()
229
230 call this%x_dof%free()
231 call this%y_dof%free()
232 call this%z_dof%free()
233
234 call this%x_interface_dof%free()
235 call this%y_interface_dof%free()
236 call this%z_interface_dof%free()
237 call this%s_interface%free()
238
239 if (allocated(this%field_name)) then
240 deallocate(this%field_name)
241 end if
242
243 call this%interface_interpolator%free()
244
245 call this%interface_dof_mask%free()
246 call this%domain_element_mask%free()
247
248 call this%free_base()
249 end subroutine overset_interface_free
250
255 subroutine overset_interface_apply_scalar(this, x, n, time, strong)
256 class(overset_interface_t), intent(inout) :: this
257 integer, intent(in) :: n
258 real(kind=rp), intent(inout), dimension(n) :: x
259 type(time_state_t), intent(in), optional :: time
260 logical, intent(in), optional :: strong
261 logical :: strong_
262
263 if (present(strong)) then
264 strong_ = strong
265 else
266 strong_ = .true.
267 end if
268
269 if (strong_) then
270 if (.not. this%updated) then
271 call this%update(time)
272 this%updated = .true.
273 end if
274
275 call masked_copy_0(x, this%bc_s%field_bc%x, this%msk, n, this%msk(0))
276 end if
277
278 end subroutine overset_interface_apply_scalar
279
284 subroutine overset_interface_apply_scalar_dev(this, x_d, time, strong, strm)
285 class(overset_interface_t), intent(inout), target :: this
286 type(c_ptr), intent(inout) :: x_d
287 type(time_state_t), intent(in), optional :: time
288 logical, intent(in), optional :: strong
289 type(c_ptr), intent(inout) :: strm
290 logical :: strong_
291
292 if (present(strong)) then
293 strong_ = strong
294 else
295 strong_ = .true.
296 end if
297
298 if (strong_) then
299 if (.not. this%updated) then
300 call this%update(time)
301 this%updated = .true.
302 end if
303
304 if (this%msk(0) .gt. 0) then
305 call device_masked_copy_0(x_d, this%bc_s%field_bc%x_d, this%bc_s%msk_d, &
306 this%bc_s%dof%size(), this%msk(0), strm)
307 end if
308 end if
309
311
313 subroutine overset_interface_apply_vector(this, x, y, z, n, time, strong)
314 class(overset_interface_t), intent(inout) :: this
315 integer, intent(in) :: n
316 real(kind=rp), intent(inout), dimension(n) :: x
317 real(kind=rp), intent(inout), dimension(n) :: y
318 real(kind=rp), intent(inout), dimension(n) :: z
319 type(time_state_t), intent(in), optional :: time
320 logical, intent(in), optional :: strong
321
322 call neko_error("overset_interface cannot apply vector BCs.&
323 & Use overset_interface_vector instead!")
324
325 end subroutine overset_interface_apply_vector
326
328 subroutine overset_interface_apply_vector_dev(this, x_d, y_d, z_d, time, &
329 strong, strm)
330 class(overset_interface_t), intent(inout), target :: this
331 type(c_ptr), intent(inout) :: x_d
332 type(c_ptr), intent(inout) :: y_d
333 type(c_ptr), intent(inout) :: z_d
334 type(time_state_t), intent(in), optional :: time
335 logical, intent(in), optional :: strong
336 type(c_ptr), intent(inout) :: strm
337
338 call neko_error("overset_interface cannot apply vector BCs.&
339 & Use overset_interface_vector instead!")
340
342
344 subroutine overset_interface_finalize(this, only_facets)
345 class(overset_interface_t), target, intent(inout) :: this
346 logical, optional, intent(in) :: only_facets
347 logical :: only_facets_
348
349 if (present(only_facets)) then
350 only_facets_ = only_facets
351 else
352 only_facets_ = .false.
353 end if
354
355 call this%finalize_base(only_facets_)
356
357 call this%bc_s%mark_facets(this%marked_facet)
358 call this%bc_s%finalize(only_facets_)
359
360 call this%build_masks_()
361
362 call this%x_interface_dof%init(this%interface_dof_mask%size(), 'x_interface')
363 call this%y_interface_dof%init(this%interface_dof_mask%size(), 'y_interface')
364 call this%z_interface_dof%init(this%interface_dof_mask%size(), 'z_interface')
365 call this%gather_interface_dofs_()
366
367 call this%setup_interpolator_()
368
369 call this%s_interface%init(this%interface_dof_mask%size(), 's_interface')
370
371 call this%interface_dof%init(3)
372 call this%interface_dof%assign_to_vector(1, this%x_interface_dof)
373 call this%interface_dof%assign_to_vector(2, this%y_interface_dof)
374 call this%interface_dof%assign_to_vector(3, this%z_interface_dof)
375
376 call this%interface_field%init(1)
377 call this%interface_field%assign_to_vector(1, this%s_interface)
378
379 end subroutine overset_interface_finalize
380
382 subroutine overset_interface_update(this, time)
383 class(overset_interface_t), intent(inout) :: this
384 type(time_state_t), intent(in) :: time
385 type(field_t), pointer :: s
386
388 call this%morph_interface(this%interface_dof, this%interface_field, &
389 this%interface_dof_mask, time, this%name, &
390 this%find_interface)
391
393 if (this%find_interface) then
394
395 ! sync
396 call this%x_interface_dof%copy_from(device_to_host, sync = .false.)
397 call this%y_interface_dof%copy_from(device_to_host, sync = .false.)
398 call this%z_interface_dof%copy_from(device_to_host, sync = .true.)
399
400 call this%interface_interpolator%find_points(this%x_interface_dof%x, &
401 this%y_interface_dof%x, this%z_interface_dof%x, &
402 this%x_interface_dof%size())
403 this%find_interface = .false.
404
405 end if
406
407 s => neko_registry%get_field(trim(this%field_name))
408
409 call this%interface_interpolator%evaluate_masked(this%s_interface%x, s%x, &
410 this%domain_element_mask, .false.)
411
412 call vector_masked_scatter_copy(this%bc_s%field_bc%x(:,1,1,1), this%s_interface, &
413 this%interface_dof_mask, this%bc_s%dof%size())
414
415 end subroutine overset_interface_update
416
417 !===================
418 ! Helper subroutines
419 !===================
420
422 subroutine build_masks_(this)
423 class(overset_interface_t), intent(inout) :: this
424 type(mask_t) :: temp_mask
425 logical, allocatable :: found(:)
426 integer :: i, j, k, e, nelems
427 integer :: lx, ly, lz
428 integer :: nonlinear_idx(4), linear_idx
429 type(stack_i4_t) :: idx_stack
430
431 call this%interface_dof_mask%init(this%msk(1:this%msk(0)), this%msk(0))
432
433 lx = this%Xh%lx
434 ly = this%Xh%ly
435 lz = this%Xh%lz
436
437 allocate(found(this%msh%nelv))
438 found = .false.
439
440 do i = 1, this%msk(0)
441 linear_idx = this%msk(i)
442 nonlinear_idx = nonlinear_index(linear_idx, lx, ly, lz)
443 found(nonlinear_idx(4)) = .true.
444 end do
445
446 nelems = 0
447 call idx_stack%init()
448 do e = 1, this%msh%nelv
449 if (found(e)) then
450 nelems = nelems + 1
451 do k = 1, this%Xh%lz
452 do j = 1, this%Xh%ly
453 do i = 1, this%Xh%lx
454 linear_idx = linear_index(i, j, k, e, lx, ly, lz)
455 call idx_stack%push(linear_idx)
456 end do
457 end do
458 end do
459 end if
460 end do
461
462 deallocate(found)
463
464 call temp_mask%init(idx_stack%array(), idx_stack%size())
465 call idx_stack%free()
466
467 call this%domain_element_mask%invert_mask(temp_mask, this%dof%size())
468 call temp_mask%free()
469
470 end subroutine build_masks_
471
473 subroutine gather_interface_dofs_(this)
474 class(overset_interface_t), intent(inout) :: this
475
476 call vector_masked_gather_copy(this%x_interface_dof, this%dof%x(:,1,1,1), &
477 this%interface_dof_mask, this%dof%size())
478 call vector_masked_gather_copy(this%y_interface_dof, this%dof%y(:,1,1,1), &
479 this%interface_dof_mask, this%dof%size())
480 call vector_masked_gather_copy(this%z_interface_dof, this%dof%z(:,1,1,1), &
481 this%interface_dof_mask, this%dof%size())
482
483 call this%x_interface_dof%copy_from(device_to_host, sync = .false.)
484 call this%y_interface_dof%copy_from(device_to_host, sync = .false.)
485 call this%z_interface_dof%copy_from(device_to_host, sync = .true.)
486
487 end subroutine gather_interface_dofs_
488
490 subroutine setup_interpolator_(this)
491 class(overset_interface_t), intent(inout) :: this
492
493 call this%interface_interpolator%init(this%dof, &
495 tol = this%interpolation_settings%tolerance, &
496 pad = this%interpolation_settings%padding, &
497 mask = this%domain_element_mask)
498
499 call this%interface_interpolator%find_points(this%x_interface_dof%x, &
500 this%y_interface_dof%x, this%z_interface_dof%x, &
501 this%x_interface_dof%size())
502
503 end subroutine setup_interpolator_
504
505end module overset_interface
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
User callback for overset-interface morphing and boundary-value updates.
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 device_to_host
Definition device.F90:47
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.
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
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 scalar boundary conditions.
subroutine overset_interface_finalize(this, only_facets)
Finalize by building the mask arrays and preparing interpolation data.
subroutine overset_interface_update(this, time)
Update values at the overset interface.
subroutine overset_interface_apply_scalar(this, x, n, time, strong)
Apply scalar.
subroutine gather_interface_dofs_(this)
Gather interface dofs.
subroutine overset_interface_apply_vector(this, x, y, z, n, time, strong)
(No-op) Apply vector.
subroutine overset_interface_free(this)
Destructor.
subroutine overset_interface_apply_scalar_dev(this, x_d, time, strong, strm)
Apply scalar (device).
subroutine overset_interface_init(this, coef, json)
Constructor.
subroutine overset_interface_init_from_components(this, coef, field_name, tol, pad)
Constructor from components.
subroutine setup_interpolator_(this)
Set up the global interpolator.
subroutine build_masks_(this)
Build masks.
subroutine overset_interface_apply_vector_dev(this, x_d, y_d, z_d, time, strong, strm)
(No-op) Apply vector (device).
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
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_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
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
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
Overset interface BC for a scalar field.
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