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
49 use vector_list, only : vector_list_t
52 use device, only : device_to_host
56 use stack, only : stack_i4_t
57 use json_module, only : json_file
59 use field, only : field_t
60 use logger, only : neko_log, log_size
62 use mpi_f08, only : mpi_allreduce, mpi_integer, mpi_sum
63 use, intrinsic :: iso_c_binding, only : c_ptr
64 use time_state, only : time_state_t
65 implicit none
66 private
67
69 type, public, extends(bc_t) :: overset_interface_t
71 type(field_dirichlet_t) :: bc_s
75 character(len=:), allocatable :: field_name
77 type(global_interpolation_t) :: interface_interpolator
79 type(mask_t) :: interface_dof_mask
80 type(mask_t) :: domain_element_mask
82 type(vector_t) :: x_dof, y_dof, z_dof
83 type(vector_t) :: x_interface_dof, y_interface_dof, z_interface_dof
85 type(vector_t) :: s_interface
86 type(vector_series_t) :: s_interface_lag
87 integer :: iextm_order = 1
88 integer :: last_tstep = -1
89 type(vector_list_t) :: interface_dof, interface_field
91 type(global_interpolation_settings_t) :: interpolation_settings
92 integer :: n_int_tot = 0
93 logical :: find_interface = .false.
94 logical :: setup = .false.
95 logical :: log = .false.
96
99 procedure(morph_overset_interface), nopass, pointer :: morph_interface => null()
100
101 contains
103 procedure, pass(this) :: init => overset_interface_init
105 procedure, pass(this) :: init_from_components => &
108 procedure, pass(this) :: free => overset_interface_free
110 procedure, pass(this) :: finalize => overset_interface_finalize
112 procedure, pass(this) :: apply_scalar => overset_interface_apply_scalar
114 procedure, pass(this) :: apply_vector => overset_interface_apply_vector
116 procedure, pass(this) :: apply_vector_dev => overset_interface_apply_vector_dev
118 procedure, pass(this) :: apply_scalar_dev => overset_interface_apply_scalar_dev
119 procedure, pass(this) :: update => overset_interface_update
120
122 procedure, pass(this), private :: build_masks_ => build_masks_
124 procedure, pass(this), private :: gather_interface_dofs_ => gather_interface_dofs_
126 procedure, pass(this), private :: setup_interpolator_ => setup_interpolator_
128 procedure, pass(this), private :: log_interface_error_ => log_interface_error_
129 end type overset_interface_t
130
131 abstract interface
132
147 subroutine morph_overset_interface(interface_dof, interface_field, &
148 interface_mask, time, bc_name, &
149 find_interface)
151 type(vector_list_t), intent(inout) :: interface_dof
152 type(vector_list_t), intent(inout) :: interface_field
153 type(mask_t), intent(in) :: interface_mask
154 type(time_state_t), intent(in) :: time
155 character(len=*), intent(in) :: bc_name
156 logical, intent(inout) :: find_interface
157 end subroutine morph_overset_interface
158 end interface
159
161
162contains
163
167 subroutine overset_interface_init(this, coef, json)
168 class(overset_interface_t), intent(inout), target :: this
169 type(coef_t), target, intent(in) :: coef
170 type(json_file), intent(inout) :: json
171 character(len=:), allocatable :: field_name
172 real(kind=rp) :: tol, pad
173 logical :: log
174
175 call json_get(json, "field_name", field_name)
176 call json_get_or_default(json, "interpolation.tolerance", tol, -1.0_rp)
177 call json_get_or_default(json, "interpolation.padding", pad, -1.0_rp)
178 call json_get_or_default(json, "order", this%iextm_order, 1)
179 if (this%iextm_order .lt. 1 .or. this%iextm_order .gt. 3) then
180 call neko_error("The order of the IEXTm time scheme must be 1 to 3.")
181 end if
182 call json_get_or_default(json, "log", log, .false.)
183
184 call this%init_from_components(coef, field_name, tol, pad, log)
185 if (allocated(field_name)) deallocate(field_name)
186
187 end subroutine overset_interface_init
188
191 subroutine overset_interface_init_from_components(this, coef, field_name, tol, pad, log)
192 class(overset_interface_t), intent(inout), target :: this
193 type(coef_t), intent(in) :: coef
194 character(len=*), intent(in) :: field_name
195 real(kind=rp), intent(in), optional :: tol, pad
196 logical, intent(in), optional :: log
197 character(len=256) :: log_buf
198
199 call this%init_base(coef)
200
201 if (present(tol)) then
202 if (tol .gt. 0.0_rp) then
203 this%interpolation_settings%tolerance = tol
204 end if
205 end if
206
207 if (present(pad)) then
208 if (pad .gt. 0.0_rp) then
209 this%interpolation_settings%padding = pad
210 end if
211 end if
212
213 if (present(log)) then
214 this%log = log
215 end if
216
217 this%field_name = field_name
218 write (log_buf, '(A,A)') "Coupling overset interface for: ", trim(this%field_name)
219 call neko_log%message(log_buf)
220
221 call this%bc_s%init_from_components(coef, this%field_name)
222 call this%field_list%init(1)
223 call this%field_list%assign_to_field(1, this%bc_s%field_bc)
224
225 call this%x_dof%init(this%dof%size(), 'x')
226 call this%y_dof%init(this%dof%size(), 'y')
227 call this%z_dof%init(this%dof%size(), 'z')
228
229 if (neko_bcknd_device .eq. 1) then
230 call device_copy(this%x_dof%x_d, this%dof%x_d, this%dof%size())
231 call device_copy(this%y_dof%x_d, this%dof%y_d, this%dof%size())
232 call device_copy(this%z_dof%x_d, this%dof%z_d, this%dof%size())
233
234 call this%x_dof%copy_from(device_to_host, sync = .false.)
235 call this%y_dof%copy_from(device_to_host, sync = .false.)
236 call this%z_dof%copy_from(device_to_host, sync = .true.)
237 else
238 call copy(this%x_dof%x, this%dof%x, this%dof%size())
239 call copy(this%y_dof%x, this%dof%y, this%dof%size())
240 call copy(this%z_dof%x, this%dof%z, this%dof%size())
241 end if
242
244
246 subroutine overset_interface_free(this)
247 class(overset_interface_t), target, intent(inout) :: this
248
249 call this%bc_s%free()
250 call this%field_list%free()
251 call this%interface_dof%free()
252 call this%interface_field%free()
253
254 call this%x_dof%free()
255 call this%y_dof%free()
256 call this%z_dof%free()
257
258 call this%x_interface_dof%free()
259 call this%y_interface_dof%free()
260 call this%z_interface_dof%free()
261 call this%s_interface%free()
262 call this%s_interface_lag%free()
263
264 if (allocated(this%field_name)) then
265 deallocate(this%field_name)
266 end if
267
268 call this%interface_interpolator%free()
269
270 call this%interface_dof_mask%free()
271 call this%domain_element_mask%free()
272
273 call this%free_base()
274 end subroutine overset_interface_free
275
280 subroutine overset_interface_apply_scalar(this, x, n, time, strong)
281 class(overset_interface_t), intent(inout) :: this
282 integer, intent(in) :: n
283 real(kind=rp), intent(inout), dimension(n) :: x
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 if (.not. this%updated) then
296 call this%update(time)
297 this%updated = .true.
298 end if
299
300 call masked_copy_0(x, this%bc_s%field_bc%x, this%msk, n, this%msk(0))
301 end if
302
303 end subroutine overset_interface_apply_scalar
304
309 subroutine overset_interface_apply_scalar_dev(this, x_d, time, strong, strm)
310 class(overset_interface_t), intent(inout), target :: this
311 type(c_ptr), intent(inout) :: x_d
312 type(time_state_t), intent(in), optional :: time
313 logical, intent(in), optional :: strong
314 type(c_ptr), intent(inout) :: strm
315 logical :: strong_
316
317 if (present(strong)) then
318 strong_ = strong
319 else
320 strong_ = .true.
321 end if
322
323 if (strong_) then
324 if (.not. this%updated) then
325 call this%update(time)
326 this%updated = .true.
327 end if
328
329 if (this%msk(0) .gt. 0) then
330 call device_masked_copy_0(x_d, this%bc_s%field_bc%x_d, this%bc_s%msk_d, &
331 this%bc_s%dof%size(), this%msk(0), strm)
332 end if
333 end if
334
336
338 subroutine overset_interface_apply_vector(this, x, y, z, n, time, strong)
339 class(overset_interface_t), intent(inout) :: this
340 integer, intent(in) :: n
341 real(kind=rp), intent(inout), dimension(n) :: x
342 real(kind=rp), intent(inout), dimension(n) :: y
343 real(kind=rp), intent(inout), dimension(n) :: z
344 type(time_state_t), intent(in), optional :: time
345 logical, intent(in), optional :: strong
346
347 call neko_error("overset_interface cannot apply vector BCs.&
348 & Use overset_interface_vector instead!")
349
350 end subroutine overset_interface_apply_vector
351
353 subroutine overset_interface_apply_vector_dev(this, x_d, y_d, z_d, time, &
354 strong, strm)
355 class(overset_interface_t), intent(inout), target :: this
356 type(c_ptr), intent(inout) :: x_d
357 type(c_ptr), intent(inout) :: y_d
358 type(c_ptr), intent(inout) :: z_d
359 type(time_state_t), intent(in), optional :: time
360 logical, intent(in), optional :: strong
361 type(c_ptr), intent(inout) :: strm
362
363 call neko_error("overset_interface cannot apply vector BCs.&
364 & Use overset_interface_vector instead!")
365
367
369 subroutine overset_interface_finalize(this, only_facets)
370 class(overset_interface_t), target, intent(inout) :: this
371 logical, optional, intent(in) :: only_facets
372 logical :: only_facets_
373
374 if (present(only_facets)) then
375 only_facets_ = only_facets
376 else
377 only_facets_ = .false.
378 end if
379
380 call this%finalize_base(only_facets_)
381
382 call this%bc_s%mark_facets(this%marked_facet)
383 call this%bc_s%finalize(only_facets_)
384
385 call this%build_masks_()
386
387 call this%x_interface_dof%init(this%interface_dof_mask%size(), 'x_interface')
388 call this%y_interface_dof%init(this%interface_dof_mask%size(), 'y_interface')
389 call this%z_interface_dof%init(this%interface_dof_mask%size(), 'z_interface')
390 call this%gather_interface_dofs_()
391
392 call this%setup_interpolator_()
393
394 call this%s_interface%init(this%interface_dof_mask%size(), 's_interface')
395
396 call this%interface_dof%init(3)
397 call this%interface_dof%assign_to_vector(1, this%x_interface_dof)
398 call this%interface_dof%assign_to_vector(2, this%y_interface_dof)
399 call this%interface_dof%assign_to_vector(3, this%z_interface_dof)
400
401 call this%interface_field%init(1)
402 call this%interface_field%assign_to_vector(1, this%s_interface)
403
404 call this%s_interface_lag%init(this%s_interface, this%iextm_order)
405
406 call mpi_allreduce(this%s_interface%size(), this%n_int_tot, 1, mpi_integer, &
407 mpi_sum, neko_global_comm)
408
409 end subroutine overset_interface_finalize
410
412 subroutine overset_interface_update(this, time)
413 class(overset_interface_t), intent(inout) :: this
414 type(time_state_t), intent(in) :: time
415 type(field_t), pointer :: s
416 type(iextm_time_scheme_t) :: time_scheme
417 integer :: nhist, ihist
418 real(kind=rp) :: iextm_coeffs(4)
419
421 call this%morph_interface(this%interface_dof, this%interface_field, &
422 this%interface_dof_mask, time, this%name, &
423 this%find_interface)
424
426 if (this%find_interface) then
427
428 ! sync
429 call this%x_interface_dof%copy_from(device_to_host, sync = .false.)
430 call this%y_interface_dof%copy_from(device_to_host, sync = .false.)
431 call this%z_interface_dof%copy_from(device_to_host, sync = .true.)
432
433 call this%interface_interpolator%find_points(this%x_interface_dof%x, &
434 this%y_interface_dof%x, this%z_interface_dof%x, &
435 this%x_interface_dof%size())
436 this%find_interface = .false.
437
438 end if
439
440 s => neko_registry%get_field(trim(this%field_name))
441
442 call this%interface_interpolator%evaluate_masked(this%s_interface%x, s%x, &
443 this%domain_element_mask, .false.)
444
445 if (this%log) then
446 call this%log_interface_error_(s)
447 end if
448
449 if (time%tstep .ne. this%last_tstep) then
450 this%last_tstep = time%tstep
451
452 call this%s_interface_lag%update()
453
454 nhist = min(time%tstep, this%iextm_order)
455 call time_scheme%compute_coeffs(iextm_coeffs, time%dtlag, nhist)
456
457 call vector_cmult2(this%s_interface, this%s_interface_lag%lv(1), iextm_coeffs(1))
458 do ihist = 2, nhist
459 call vector_add2s2(this%s_interface, this%s_interface_lag%lv(ihist), iextm_coeffs(ihist))
460 end do
461 end if
462
463 call vector_masked_scatter_copy(this%bc_s%field_bc%x(:,1,1,1), this%s_interface, &
464 this%interface_dof_mask, this%bc_s%dof%size())
465
466 end subroutine overset_interface_update
467
469 subroutine log_interface_error_(this, s)
470 class(overset_interface_t), intent(inout) :: this
471 type(field_t), pointer, intent(in) :: s
472 real(kind=rp) :: s_int_norm
473 type(vector_t), pointer :: error
474 integer :: ind(1)
475 logical :: clear_scratch = .false.
476 character(len=256) :: log_buf
477
478 call neko_scratch_registry%request_vector(error, ind(1), this%s_interface%size(), &
479 clear_scratch)
480 call vector_masked_gather_copy(error, s%x(:,1,1,1), this%interface_dof_mask, &
481 this%dof%size())
482 call vector_add2s2(error, this%s_interface, -1.0_rp)
483 s_int_norm = sqrt(vector_glsc2(error, error)) / sqrt(real(this%n_int_tot, kind=rp))
484 call neko_scratch_registry%relinquish(ind)
485
486 write(log_buf, '(A12,A3,A10,1x,E15.7)') 'Interface BC', ' | ', &
487 'L2 Error: ', s_int_norm
488 call neko_log%message(log_buf)
489
490 end subroutine log_interface_error_
491
492 !===================
493 ! Helper subroutines
494 !===================
495
497 subroutine build_masks_(this)
498 class(overset_interface_t), intent(inout) :: this
499 type(mask_t) :: temp_mask
500 logical, allocatable :: found(:)
501 integer :: i, j, k, e, nelems
502 integer :: lx, ly, lz
503 integer :: nonlinear_idx(4), linear_idx
504 type(stack_i4_t) :: idx_stack
505
506 call this%interface_dof_mask%init(this%msk(1:this%msk(0)), this%msk(0))
507
508 lx = this%Xh%lx
509 ly = this%Xh%ly
510 lz = this%Xh%lz
511
512 allocate(found(this%msh%nelv))
513 found = .false.
514
515 do i = 1, this%msk(0)
516 linear_idx = this%msk(i)
517 nonlinear_idx = nonlinear_index(linear_idx, lx, ly, lz)
518 found(nonlinear_idx(4)) = .true.
519 end do
520
521 nelems = 0
522 call idx_stack%init()
523 do e = 1, this%msh%nelv
524 if (found(e)) then
525 nelems = nelems + 1
526 do k = 1, this%Xh%lz
527 do j = 1, this%Xh%ly
528 do i = 1, this%Xh%lx
529 linear_idx = linear_index(i, j, k, e, lx, ly, lz)
530 call idx_stack%push(linear_idx)
531 end do
532 end do
533 end do
534 end if
535 end do
536
537 deallocate(found)
538
539 call temp_mask%init(idx_stack%array(), idx_stack%size())
540 call idx_stack%free()
541
542 call this%domain_element_mask%invert_mask(temp_mask, this%dof%size())
543 call temp_mask%free()
544
545 end subroutine build_masks_
546
548 subroutine gather_interface_dofs_(this)
549 class(overset_interface_t), intent(inout) :: this
550
551 call vector_masked_gather_copy(this%x_interface_dof, this%dof%x(:,1,1,1), &
552 this%interface_dof_mask, this%dof%size())
553 call vector_masked_gather_copy(this%y_interface_dof, this%dof%y(:,1,1,1), &
554 this%interface_dof_mask, this%dof%size())
555 call vector_masked_gather_copy(this%z_interface_dof, this%dof%z(:,1,1,1), &
556 this%interface_dof_mask, this%dof%size())
557
558 call this%x_interface_dof%copy_from(device_to_host, sync = .false.)
559 call this%y_interface_dof%copy_from(device_to_host, sync = .false.)
560 call this%z_interface_dof%copy_from(device_to_host, sync = .true.)
561
562 end subroutine gather_interface_dofs_
563
565 subroutine setup_interpolator_(this)
566 class(overset_interface_t), intent(inout) :: this
567
568 call this%interface_interpolator%init(this%dof, &
570 tol = this%interpolation_settings%tolerance, &
571 pad = this%interpolation_settings%padding, &
572 mask = this%domain_element_mask)
573
574 call this%interface_interpolator%find_points(this%x_interface_dof%x, &
575 this%y_interface_dof%x, this%z_interface_dof%x, &
576 this%x_interface_dof%size())
577
578 end subroutine setup_interpolator_
579
580end module overset_interface
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
double real
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:46
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:48
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:80
integer, parameter, public log_size
Definition log.f90:46
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:289
subroutine, public masked_copy_0(a, b, mask, n, n_mask)
Copy a masked vector .
Definition math.f90:310
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 log_interface_error_(this, s)
Log interface RMSE for the scalar field.
subroutine overset_interface_init(this, coef, json)
Constructor.
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).
subroutine overset_interface_init_from_components(this, coef, field_name, tol, pad, log)
Constructor from components.
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
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Implements a dynamic stack ADT.
Definition stack.f90:49
Base class for time integration schemes.
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:289
subroutine, public vector_masked_gather_copy(a, b, mask, n)
Gather a vector to reduced contigous array .
real(kind=rp) function, public vector_glsc2(a, b, n)
subroutine, public vector_add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
subroutine, public vector_masked_scatter_copy(a, b, mask, n)
Scatter a contiguous vector into an array .
subroutine, public vector_cmult2(a, b, c, n)
Multiplication by constant c .
Contains the vector_series_t type.
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.
Explicit interface extrapolation scheme for overset grids.
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
Stores a series (sequence) of vectors, logically connected to a base vector, and arranged according t...