54 use json_module,
only : json_file
58 use,
intrinsic :: iso_c_binding, only : c_ptr
70 character(len=:),
allocatable :: field_name
78 type(
vector_t) :: x_interface_dof, y_interface_dof, z_interface_dof
84 logical :: find_interface = .false.
85 logical :: setup = .false.
95 procedure, pass(this) :: init_from_components => &
136 interface_mask, time, bc_name, &
141 type(
mask_t),
intent(in) :: interface_mask
143 character(len=*),
intent(in) :: bc_name
144 logical,
intent(inout) :: find_interface
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
162 call json_get(json,
"field_name", field_name)
165 call this%init_from_components(coef, field_name, tol, pad)
166 if (
allocated(field_name))
deallocate(field_name)
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
179 call this%init_base(coef)
181 if (
present(tol))
then
182 if (tol .gt. 0.0_rp)
then
183 this%interpolation_settings%tolerance = tol
187 if (
present(pad))
then
188 if (pad .gt. 0.0_rp)
then
189 this%interpolation_settings%padding = pad
193 this%field_name = field_name
194 write (log_buf,
'(A,A)')
"Coupling overset interface for: ", trim(this%field_name)
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)
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')
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())
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())
225 call this%bc_s%free()
226 call this%field_list%free()
227 call this%interface_dof%free()
228 call this%interface_field%free()
230 call this%x_dof%free()
231 call this%y_dof%free()
232 call this%z_dof%free()
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()
239 if (
allocated(this%field_name))
then
240 deallocate(this%field_name)
243 call this%interface_interpolator%free()
245 call this%interface_dof_mask%free()
246 call this%domain_element_mask%free()
248 call this%free_base()
257 integer,
intent(in) :: n
258 real(kind=
rp),
intent(inout),
dimension(n) :: x
260 logical,
intent(in),
optional :: strong
263 if (
present(strong))
then
270 if (.not. this%updated)
then
271 call this%update(time)
272 this%updated = .true.
275 call masked_copy_0(x, this%bc_s%field_bc%x, this%msk, n, this%msk(0))
286 type(c_ptr),
intent(inout) :: x_d
288 logical,
intent(in),
optional :: strong
289 type(c_ptr),
intent(inout) :: strm
292 if (
present(strong))
then
299 if (.not. this%updated)
then
300 call this%update(time)
301 this%updated = .true.
304 if (this%msk(0) .gt. 0)
then
306 this%bc_s%dof%size(), this%msk(0), strm)
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
320 logical,
intent(in),
optional :: strong
322 call neko_error(
"overset_interface cannot apply vector BCs.&
323 & Use overset_interface_vector instead!")
331 type(c_ptr),
intent(inout) :: x_d
332 type(c_ptr),
intent(inout) :: y_d
333 type(c_ptr),
intent(inout) :: z_d
335 logical,
intent(in),
optional :: strong
336 type(c_ptr),
intent(inout) :: strm
338 call neko_error(
"overset_interface cannot apply vector BCs.&
339 & Use overset_interface_vector instead!")
346 logical,
optional,
intent(in) :: only_facets
347 logical :: only_facets_
349 if (
present(only_facets))
then
350 only_facets_ = only_facets
352 only_facets_ = .false.
355 call this%finalize_base(only_facets_)
357 call this%bc_s%mark_facets(this%marked_facet)
358 call this%bc_s%finalize(only_facets_)
360 call this%build_masks_()
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_()
367 call this%setup_interpolator_()
369 call this%s_interface%init(this%interface_dof_mask%size(),
's_interface')
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)
376 call this%interface_field%init(1)
377 call this%interface_field%assign_to_vector(1, this%s_interface)
388 call this%morph_interface(this%interface_dof, this%interface_field, &
389 this%interface_dof_mask, time, this%name, &
393 if (this%find_interface)
then
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.)
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.
409 call this%interface_interpolator%evaluate_masked(this%s_interface%x, s%x, &
410 this%domain_element_mask, .false.)
413 this%interface_dof_mask, this%bc_s%dof%size())
425 logical,
allocatable :: found(:)
426 integer :: i, j, k, e, nelems
427 integer :: lx, ly, lz
428 integer :: nonlinear_idx(4), linear_idx
431 call this%interface_dof_mask%init(this%msk(1:this%msk(0)), this%msk(0))
437 allocate(found(this%msh%nelv))
440 do i = 1, this%msk(0)
441 linear_idx = this%msk(i)
443 found(nonlinear_idx(4)) = .true.
447 call idx_stack%init()
448 do e = 1, this%msh%nelv
455 call idx_stack%push(linear_idx)
464 call temp_mask%init(idx_stack%array(), idx_stack%size())
465 call idx_stack%free()
467 call this%domain_element_mask%invert_mask(temp_mask, this%dof%size())
468 call temp_mask%free()
477 this%interface_dof_mask, this%dof%size())
479 this%interface_dof_mask, this%dof%size())
481 this%interface_dof_mask, this%dof%size())
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.)
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)
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())
__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 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.
type(registry_t), target, public neko_registry
Global field registry.
Implements a dynamic stack ADT.
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 .
subroutine, public vector_masked_scatter_copy(a, b, mask, n)
Scatter a contiguous vector into an array .
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.
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