57 use json_module,
only : json_file
62 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_sum
63 use,
intrinsic :: iso_c_binding, only : c_ptr
75 character(len=:),
allocatable :: field_name
83 type(
vector_t) :: x_interface_dof, y_interface_dof, z_interface_dof
87 integer :: iextm_order = 1
88 integer :: last_tstep = -1
92 integer :: n_int_tot = 0
93 logical :: find_interface = .false.
94 logical :: setup = .false.
95 logical :: log = .false.
105 procedure, pass(this) :: init_from_components => &
148 interface_mask, time, bc_name, &
153 type(
mask_t),
intent(in) :: interface_mask
155 character(len=*),
intent(in) :: bc_name
156 logical,
intent(inout) :: find_interface
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
175 call json_get(json,
"field_name", field_name)
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.")
184 call this%init_from_components(coef, field_name, tol, pad, log)
185 if (
allocated(field_name))
deallocate(field_name)
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
199 call this%init_base(coef)
201 if (
present(tol))
then
202 if (tol .gt. 0.0_rp)
then
203 this%interpolation_settings%tolerance = tol
207 if (
present(pad))
then
208 if (pad .gt. 0.0_rp)
then
209 this%interpolation_settings%padding = pad
213 if (
present(log))
then
217 this%field_name = field_name
218 write (log_buf,
'(A,A)')
"Coupling overset interface for: ", trim(this%field_name)
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)
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')
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())
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())
249 call this%bc_s%free()
250 call this%field_list%free()
251 call this%interface_dof%free()
252 call this%interface_field%free()
254 call this%x_dof%free()
255 call this%y_dof%free()
256 call this%z_dof%free()
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()
264 if (
allocated(this%field_name))
then
265 deallocate(this%field_name)
268 call this%interface_interpolator%free()
270 call this%interface_dof_mask%free()
271 call this%domain_element_mask%free()
273 call this%free_base()
282 integer,
intent(in) :: n
283 real(kind=
rp),
intent(inout),
dimension(n) :: x
285 logical,
intent(in),
optional :: strong
288 if (
present(strong))
then
295 if (.not. this%updated)
then
296 call this%update(time)
297 this%updated = .true.
300 call masked_copy_0(x, this%bc_s%field_bc%x, this%msk, n, this%msk(0))
311 type(c_ptr),
intent(inout) :: x_d
313 logical,
intent(in),
optional :: strong
314 type(c_ptr),
intent(inout) :: strm
317 if (
present(strong))
then
324 if (.not. this%updated)
then
325 call this%update(time)
326 this%updated = .true.
329 if (this%msk(0) .gt. 0)
then
331 this%bc_s%dof%size(), this%msk(0), strm)
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
345 logical,
intent(in),
optional :: strong
347 call neko_error(
"overset_interface cannot apply vector BCs.&
348 & Use overset_interface_vector instead!")
356 type(c_ptr),
intent(inout) :: x_d
357 type(c_ptr),
intent(inout) :: y_d
358 type(c_ptr),
intent(inout) :: z_d
360 logical,
intent(in),
optional :: strong
361 type(c_ptr),
intent(inout) :: strm
363 call neko_error(
"overset_interface cannot apply vector BCs.&
364 & Use overset_interface_vector instead!")
371 logical,
optional,
intent(in) :: only_facets
372 logical :: only_facets_
374 if (
present(only_facets))
then
375 only_facets_ = only_facets
377 only_facets_ = .false.
380 call this%finalize_base(only_facets_)
382 call this%bc_s%mark_facets(this%marked_facet)
383 call this%bc_s%finalize(only_facets_)
385 call this%build_masks_()
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_()
392 call this%setup_interpolator_()
394 call this%s_interface%init(this%interface_dof_mask%size(),
's_interface')
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)
401 call this%interface_field%init(1)
402 call this%interface_field%assign_to_vector(1, this%s_interface)
404 call this%s_interface_lag%init(this%s_interface, this%iextm_order)
406 call mpi_allreduce(this%s_interface%size(), this%n_int_tot, 1, mpi_integer, &
417 integer :: nhist, ihist
418 real(kind=
rp) :: iextm_coeffs(4)
421 call this%morph_interface(this%interface_dof, this%interface_field, &
422 this%interface_dof_mask, time, this%name, &
426 if (this%find_interface)
then
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.)
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.
442 call this%interface_interpolator%evaluate_masked(this%s_interface%x, s%x, &
443 this%domain_element_mask, .false.)
446 call this%log_interface_error_(s)
449 if (time%tstep .ne. this%last_tstep)
then
450 this%last_tstep = time%tstep
452 call this%s_interface_lag%update()
454 nhist = min(time%tstep, this%iextm_order)
455 call time_scheme%compute_coeffs(iextm_coeffs, time%dtlag, nhist)
457 call vector_cmult2(this%s_interface, this%s_interface_lag%lv(1), iextm_coeffs(1))
459 call vector_add2s2(this%s_interface, this%s_interface_lag%lv(ihist), iextm_coeffs(ihist))
464 this%interface_dof_mask, this%bc_s%dof%size())
471 type(
field_t),
pointer,
intent(in) :: s
472 real(kind=
rp) :: s_int_norm
475 logical :: clear_scratch = .false.
476 character(len=256) :: log_buf
486 write(log_buf,
'(A12,A3,A10,1x,E15.7)')
'Interface BC',
' | ', &
487 'L2 Error: ', s_int_norm
500 logical,
allocatable :: found(:)
501 integer :: i, j, k, e, nelems
502 integer :: lx, ly, lz
503 integer :: nonlinear_idx(4), linear_idx
506 call this%interface_dof_mask%init(this%msk(1:this%msk(0)), this%msk(0))
512 allocate(found(this%msh%nelv))
515 do i = 1, this%msk(0)
516 linear_idx = this%msk(i)
518 found(nonlinear_idx(4)) = .true.
522 call idx_stack%init()
523 do e = 1, this%msh%nelv
530 call idx_stack%push(linear_idx)
539 call temp_mask%init(idx_stack%array(), idx_stack%size())
540 call idx_stack%free()
542 call this%domain_element_mask%invert_mask(temp_mask, this%dof%size())
543 call temp_mask%free()
552 this%interface_dof_mask, this%dof%size())
554 this%interface_dof_mask, this%dof%size())
556 this%interface_dof_mask, this%dof%size())
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.)
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)
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())
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
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.
type(mpi_comm), public neko_global_comm
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.
integer, parameter, public device_to_host
Defines user dirichlet condition for a scalar field.
Implements global_interpolation given a dofmap.
Utilities for retrieving parameters from the case files.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
Object for handling masks in Neko.
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public masked_copy_0(a, b, mask, n, n_mask)
Copy a masked vector .
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
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.
type(registry_t), target, public neko_registry
Global field registry.
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.
Base class for time integration schemes.
Module with things related to the simulation time.
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, :)
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.
Base type for a boundary condition.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
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...
Overset interface BC for a scalar field.
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...