63 use json_module,
only : json_file
66 use,
intrinsic :: iso_c_binding, only : c_ptr, c_size_t
89 type(
vector_t) :: x_interface_dof, y_interface_dof, z_interface_dof
90 type(
vector_t) :: u_interface, v_interface, w_interface
94 logical :: find_interface = .false.
95 logical :: setup = .false.
105 procedure, pass(this) :: init_from_components => &
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 => &
141 type(
coef_t),
target,
intent(in) :: coef
142 type(json_file),
intent(inout) ::json
143 real(kind=
rp) :: tol, pad
151 call this%init_from_components(coef, tol, pad)
159 type(
coef_t),
intent(in) :: coef
160 real(kind=
rp),
intent(in),
optional :: tol, pad
163 call this%init_base(coef)
166 if (
present(tol))
then
167 if (tol .gt. 0.0_rp) this%interpolation_settings%tolerance = tol
169 if (
present(pad))
then
170 if (pad .gt. 0.0_rp) this%interpolation_settings%padding = pad
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")
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)
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')
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())
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())
211 call this%bc_u%free()
212 call this%bc_v%free()
213 call this%bc_w%free()
215 call this%field_list%free()
216 call this%interface_dof%free()
217 call this%interface_field%free()
219 call this%x_dof%free()
220 call this%y_dof%free()
221 call this%z_dof%free()
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()
245 integer,
intent(in) :: n
246 real(kind=
rp),
intent(inout),
dimension(n) :: x
248 logical,
intent(in),
optional :: strong
250 call neko_error(
"overset_interface_vector cannot apply scalar BCs.&
251 & Use overset_interface_vector::apply_vector instead!")
261 type(c_ptr),
intent(inout) :: x_d
263 logical,
intent(in),
optional :: strong
264 type(c_ptr),
intent(inout) :: strm
266 call neko_error(
"overset_interface_vector cannot apply scalar BCs.&
267 & Use overset_interface_vector::apply_vector instead!")
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
285 logical,
intent(in),
optional :: strong
288 if (
present(strong))
then
298 if (.not. this%updated)
then
299 call this%update(time)
300 this%updated = .true.
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))
320 type(c_ptr),
intent(inout) :: x_d
321 type(c_ptr),
intent(inout) :: y_d
322 type(c_ptr),
intent(inout) :: z_d
324 logical,
intent(in),
optional :: strong
325 type(c_ptr),
intent(inout) :: strm
328 if (
present(strong))
then
335 if (.not. this%updated)
then
336 call this%update(time)
337 this%updated = .true.
340 if (this%msk(0) .gt. 0)
then
342 this%bc_u%msk_d, this%bc_u%dof%size(), this%msk(0), strm)
344 this%bc_v%msk_d, this%bc_v%dof%size(), this%msk(0), strm)
346 this%bc_w%msk_d, this%bc_w%dof%size(), this%msk(0), strm)
355 logical,
optional,
intent(in) :: only_facets
356 logical :: only_facets_
358 if (
present(only_facets))
then
359 only_facets_ = only_facets
361 only_facets_ = .false.
365 call this%finalize_base(only_facets_)
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)
371 call this%bc_u%finalize(only_facets_)
372 call this%bc_v%finalize(only_facets_)
373 call this%bc_w%finalize(only_facets_)
376 call this%build_masks_()
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_()
385 call this%setup_interpolator_()
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')
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)
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)
410 type(
field_t),
pointer :: u, v, w
414 call this%morph_interface(this%interface_dof, this%interface_field, &
415 this%interface_dof_mask, time, this%name, &
425 if (this%find_interface)
then
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.)
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.
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.)
451 this%interface_dof_mask, this%bc_u%dof%size())
453 this%interface_dof_mask, this%bc_v%dof%size())
455 this%interface_dof_mask, this%bc_w%dof%size())
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
475 call this%interface_dof_mask%init(this%msk(1:this%msk(0)), this%msk(0))
482 allocate(found(this%msh%nelv))
485 do i = 1, this%msk(0)
486 linear_idx = this%msk(i)
488 found(nonlinear_idx(4)) = .true.
493 do e = 1, this%msh%nelv
500 call stack%push(linear_idx)
507 call temp_mask%init(
stack%array(),
stack%size())
511 call this%domain_element_mask%invert_mask(temp_mask, this%dof%size())
513 call temp_mask%free()
523 this%dof%x(:,1,1,1), &
524 this%interface_dof_mask, &
527 this%dof%y(:,1,1,1), &
528 this%interface_dof_mask, &
531 this%dof%z(:,1,1,1), &
532 this%interface_dof_mask, &
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.)
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)
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())
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
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 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 host_to_device
integer, parameter, public device_to_host
Defines a dirichlet boundary condition.
Defines a mapping of the degrees of freedom.
Defines user dirichlet condition for a scalar field.
Implements global_interpolation given a dofmap.
Utilities for retrieving parameters from the case files.
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 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.
type(registry_t), target, public neko_registry
Global field registry.
Implements a dynamic stack ADT.
Module with things related to the simulation time.
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,...
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_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 .
Base type for a boundary condition.
A list of allocatable `bc_t`. Follows the standard interface of lists.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Generic Dirichlet boundary condition on .
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...
Extension of the user defined dirichlet condition overset_interface
A struct that contains all info about the time, expand as needed.
vector_list_t, To be able to group vectors together