65 use json_module,
only : json_file
69 use,
intrinsic :: iso_c_binding, only : c_ptr, c_size_t
71 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_sum
97 type(
vector_t) :: x_interface_dof, y_interface_dof, z_interface_dof
98 type(
vector_t) :: u_interface, v_interface, w_interface
100 integer :: iextm_order = 1
101 integer :: last_tstep = -1
105 integer :: n_int_tot = 0
106 logical :: find_interface = .false.
107 logical :: setup = .false.
108 logical :: log = .false.
118 procedure, pass(this) :: init_from_components => &
125 procedure, pass(this) :: apply_scalar => &
128 procedure, pass(this) :: apply_vector => &
131 procedure, pass(this) :: apply_vector_dev => &
134 procedure, pass(this) :: apply_scalar_dev => &
156 type(
coef_t),
target,
intent(in) :: coef
157 type(json_file),
intent(inout) ::json
158 real(kind=
rp) :: tol, pad
167 if (this%iextm_order .lt. 1 .or. this%iextm_order .gt. 3)
then
168 call neko_error(
"The order of the IEXTm time scheme must be 1 to 3.")
172 call this%init_from_components(coef, tol, pad, log)
183 type(
coef_t),
intent(in) :: coef
184 real(kind=
rp),
intent(in),
optional :: tol, pad
185 logical,
intent(in),
optional :: log
188 call this%init_base(coef)
191 if (
present(tol))
then
192 if (tol .gt. 0.0_rp) this%interpolation_settings%tolerance = tol
194 if (
present(pad))
then
195 if (pad .gt. 0.0_rp) this%interpolation_settings%padding = pad
197 if (
present(log))
then
201 call this%bc_u%init_from_components(coef,
"u")
202 call this%bc_v%init_from_components(coef,
"v")
203 call this%bc_w%init_from_components(coef,
"w")
205 call this%field_list%init(3)
206 call this%field_list%assign_to_field(1, this%bc_u%field_bc)
207 call this%field_list%assign_to_field(2, this%bc_v%field_bc)
208 call this%field_list%assign_to_field(3, this%bc_w%field_bc)
211 call this%x_dof%init(this%dof%size(),
'x')
212 call this%y_dof%init(this%dof%size(),
'y')
213 call this%z_dof%init(this%dof%size(),
'z')
219 call device_copy(this%x_dof%x_d, this%dof%x_d, this%dof%size())
220 call device_copy(this%y_dof%x_d, this%dof%y_d, this%dof%size())
221 call device_copy(this%z_dof%x_d, this%dof%z_d, this%dof%size())
227 call copy(this%x_dof%x, this%dof%x, this%dof%size())
228 call copy(this%y_dof%x, this%dof%y, this%dof%size())
229 call copy(this%z_dof%x, this%dof%z, this%dof%size())
239 call this%bc_u%free()
240 call this%bc_v%free()
241 call this%bc_w%free()
243 call this%field_list%free()
244 call this%interface_dof%free()
245 call this%interface_field%free()
247 call this%x_dof%free()
248 call this%y_dof%free()
249 call this%z_dof%free()
251 call this%x_interface_dof%free()
252 call this%y_interface_dof%free()
253 call this%z_interface_dof%free()
254 call this%u_interface%free()
255 call this%v_interface%free()
256 call this%w_interface%free()
257 call this%u_interface_lag%free()
258 call this%v_interface_lag%free()
259 call this%w_interface_lag%free()
260 call this%interface_interpolator%free()
261 call this%interface_dof_mask%free()
262 call this%domain_element_mask%free()
263 call this%free_base()
276 integer,
intent(in) :: n
277 real(kind=
rp),
intent(inout),
dimension(n) :: x
279 logical,
intent(in),
optional :: strong
281 call neko_error(
"overset_interface_vector cannot apply scalar BCs.&
282 & Use overset_interface_vector::apply_vector instead!")
292 type(c_ptr),
intent(inout) :: x_d
294 logical,
intent(in),
optional :: strong
295 type(c_ptr),
intent(inout) :: strm
297 call neko_error(
"overset_interface_vector cannot apply scalar BCs.&
298 & Use overset_interface_vector::apply_vector instead!")
311 integer,
intent(in) :: n
312 real(kind=
rp),
intent(inout),
dimension(n) :: x
313 real(kind=
rp),
intent(inout),
dimension(n) :: y
314 real(kind=
rp),
intent(inout),
dimension(n) :: z
316 logical,
intent(in),
optional :: strong
319 if (
present(strong))
then
329 if (.not. this%updated)
then
330 call this%update(time)
331 this%updated = .true.
335 call masked_copy_0(x, this%bc_u%field_bc%x, this%msk, n, this%msk(0))
336 call masked_copy_0(y, this%bc_v%field_bc%x, this%msk, n, this%msk(0))
337 call masked_copy_0(z, this%bc_w%field_bc%x, this%msk, n, this%msk(0))
351 type(c_ptr),
intent(inout) :: x_d
352 type(c_ptr),
intent(inout) :: y_d
353 type(c_ptr),
intent(inout) :: z_d
355 logical,
intent(in),
optional :: strong
356 type(c_ptr),
intent(inout) :: strm
359 if (
present(strong))
then
366 if (.not. this%updated)
then
367 call this%update(time)
368 this%updated = .true.
371 if (this%msk(0) .gt. 0)
then
373 this%bc_u%msk_d, this%bc_u%dof%size(), this%msk(0), strm)
375 this%bc_v%msk_d, this%bc_v%dof%size(), this%msk(0), strm)
377 this%bc_w%msk_d, this%bc_w%dof%size(), this%msk(0), strm)
386 logical,
optional,
intent(in) :: only_facets
387 logical :: only_facets_
389 if (
present(only_facets))
then
390 only_facets_ = only_facets
392 only_facets_ = .false.
396 call this%finalize_base(only_facets_)
398 call this%bc_u%mark_facets(this%marked_facet)
399 call this%bc_v%mark_facets(this%marked_facet)
400 call this%bc_w%mark_facets(this%marked_facet)
402 call this%bc_u%finalize(only_facets_)
403 call this%bc_v%finalize(only_facets_)
404 call this%bc_w%finalize(only_facets_)
407 call this%build_masks_()
410 call this%x_interface_dof%init(this%interface_dof_mask%size(),
'x_interface')
411 call this%y_interface_dof%init(this%interface_dof_mask%size(),
'y_interface')
412 call this%z_interface_dof%init(this%interface_dof_mask%size(),
'z_interface')
413 call this%gather_interface_dofs_()
416 call this%setup_interpolator_()
419 call this%u_interface%init(this%interface_dof_mask%size(),
'u_interface')
420 call this%v_interface%init(this%interface_dof_mask%size(),
'v_interface')
421 call this%w_interface%init(this%interface_dof_mask%size(),
'w_interface')
424 call this%interface_dof%init(3)
425 call this%interface_dof%assign_to_vector(1, this%x_interface_dof)
426 call this%interface_dof%assign_to_vector(2, this%y_interface_dof)
427 call this%interface_dof%assign_to_vector(3, this%z_interface_dof)
429 call this%interface_field%init(3)
430 call this%interface_field%assign_to_vector(1, this%u_interface)
431 call this%interface_field%assign_to_vector(2, this%v_interface)
432 call this%interface_field%assign_to_vector(3, this%w_interface)
435 call this%u_interface_lag%init(this%u_interface, this%iextm_order)
436 call this%v_interface_lag%init(this%v_interface, this%iextm_order)
437 call this%w_interface_lag%init(this%w_interface, this%iextm_order)
439 call mpi_allreduce(this%u_interface%size(), this%n_int_tot, 1, mpi_integer, &
449 type(
field_t),
pointer :: u, v, w
451 integer :: nhist, ihist
452 real(kind=
rp) :: iextm_coeffs(4)
456 call this%morph_interface(this%interface_dof, this%interface_field, &
457 this%interface_dof_mask, time, this%name, &
467 if (this%find_interface)
then
470 call this%x_interface_dof%copy_from(
device_to_host, sync = .false.)
471 call this%y_interface_dof%copy_from(
device_to_host, sync = .false.)
472 call this%z_interface_dof%copy_from(
device_to_host, sync = .true.)
474 call this%interface_interpolator%find_points(this%x_interface_dof%x, &
475 this%y_interface_dof%x, this%z_interface_dof%x, &
476 this%x_interface_dof%size())
477 this%find_interface = .false.
487 call this%interface_interpolator%evaluate_masked(this%u_interface%x, u%x, this%domain_element_mask, .false.)
488 call this%interface_interpolator%evaluate_masked(this%v_interface%x, v%x, this%domain_element_mask, .false.)
489 call this%interface_interpolator%evaluate_masked(this%w_interface%x, w%x, this%domain_element_mask, .false.)
492 call this%log_interface_error_(u, v, w)
497 if (time%tstep .ne. this%last_tstep)
then
499 this%last_tstep = time%tstep
502 call this%u_interface_lag%update()
503 call this%v_interface_lag%update()
504 call this%w_interface_lag%update()
507 nhist = min(time%tstep, this%iextm_order)
508 call time_scheme%compute_coeffs(iextm_coeffs, time%dtlag, nhist)
511 call vector_cmult2(this%u_interface, this%u_interface_lag%lv(1), iextm_coeffs(1))
512 call vector_cmult2(this%v_interface, this%v_interface_lag%lv(1), iextm_coeffs(1))
513 call vector_cmult2(this%w_interface, this%w_interface_lag%lv(1), iextm_coeffs(1))
515 call vector_add2s2(this%u_interface, this%u_interface_lag%lv(ihist), iextm_coeffs(ihist))
516 call vector_add2s2(this%v_interface, this%v_interface_lag%lv(ihist), iextm_coeffs(ihist))
517 call vector_add2s2(this%w_interface, this%w_interface_lag%lv(ihist), iextm_coeffs(ihist))
525 this%interface_dof_mask, this%bc_u%dof%size())
527 this%interface_dof_mask, this%bc_v%dof%size())
529 this%interface_dof_mask, this%bc_w%dof%size())
537 type(
field_t),
pointer,
intent(in) :: u, v, w
538 real(kind=
rp) :: u_int_norm, v_int_norm, w_int_norm
541 logical :: clear_scratch = .false.
542 character(len=256) :: log_buf
566 write(log_buf,
'(A12,A3,A10,1x,A1,E15.7,A1,E15.7,A1,E15.7,A1)') &
567 'Interface BC',
' | ',
'L2 Error: ',
'(', &
568 u_int_norm,
',', v_int_norm,
',', w_int_norm,
')'
581 logical,
allocatable :: found(:)
582 integer :: i, j, k, e, new_size, nelems
583 integer :: lx, ly, lz
584 integer :: nonlinear_idx(4), linear_idx
588 call this%interface_dof_mask%init(this%msk(1:this%msk(0)), this%msk(0))
595 allocate(found(this%msh%nelv))
598 do i = 1, this%msk(0)
599 linear_idx = this%msk(i)
601 found(nonlinear_idx(4)) = .true.
606 do e = 1, this%msh%nelv
613 call stack%push(linear_idx)
620 call temp_mask%init(
stack%array(),
stack%size())
624 call this%domain_element_mask%invert_mask(temp_mask, this%dof%size())
626 call temp_mask%free()
636 this%dof%x(:,1,1,1), &
637 this%interface_dof_mask, &
640 this%dof%y(:,1,1,1), &
641 this%interface_dof_mask, &
644 this%dof%z(:,1,1,1), &
645 this%interface_dof_mask, &
648 call this%x_interface_dof%copy_from(
device_to_host, sync = .false.)
649 call this%y_interface_dof%copy_from(
device_to_host, sync = .false.)
650 call this%z_interface_dof%copy_from(
device_to_host, sync = .true.)
660 call this%interface_interpolator%init(this%dof, &
662 tol=this%interpolation_settings%tolerance, &
663 pad=this%interpolation_settings%padding, &
664 mask=this%domain_element_mask)
667 call this%interface_interpolator%find_points(this%x_interface_dof%x, &
668 this%y_interface_dof%x, &
669 this%z_interface_dof%x, &
670 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.
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 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_apply_scalar_dev(this, x_d, time, strong, strm)
No-op apply scalar (device).
subroutine overset_interface_vector_init_from_components(this, coef, tol, pad, log)
Constructor from components.
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 log_interface_error_(this, s)
Log interface RMSE for the scalar field.
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.
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.
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 .
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.
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.
Explicit interface extrapolation scheme for overset grids.
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
Stores a series (sequence) of vectors, logically connected to a base vector, and arranged according t...