Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
overset_interface_vector.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
40 use dirichlet, only : dirichlet_t
43 use mask, only : mask_t
44 use dofmap, only : dofmap_t
45 use bc, only : bc_t
46 use bc_list, only : bc_list_t
47 use utils, only : split_string
48 use field, only : field_t
49 use field_list, only : field_list_t
50 use math, only : masked_copy_0
52 use dofmap, only : dofmap_t
53 use vector, only : vector_t
54 use vector_list, only : vector_list_t
58 use math, only : copy
60 use vector_math, only : vector_copy
64 use stack, only: stack_i4_t
65 use json_module, only : json_file
67 use field_list, only : field_list_t
69 use, intrinsic :: iso_c_binding, only : c_ptr, c_size_t
70 use time_state, only : time_state_t
71 use mpi_f08, only : mpi_allreduce, mpi_integer, mpi_sum
73 use logger, only : neko_log, log_size
74
75
76 implicit none
77 private
78
80 ! for the application on a vector field.
81 type, public, extends(bc_t) :: overset_interface_vector_t
82 ! The bc for the first compoent.
83 type(field_dirichlet_t) :: bc_u
84 ! The bc for the second compoent.
85 type(field_dirichlet_t) :: bc_v
86 ! The bc for the third compoent.
87 type(field_dirichlet_t) :: bc_w
91 type(global_interpolation_t) :: interface_interpolator
93 type(mask_t) :: interface_dof_mask
94 type(mask_t) :: domain_element_mask
96 type(vector_t) :: x_dof, y_dof, z_dof
97 type(vector_t) :: x_interface_dof, y_interface_dof, z_interface_dof
98 type(vector_t) :: u_interface, v_interface, w_interface
99 type(vector_series_t) :: u_interface_lag, v_interface_lag, w_interface_lag
100 integer :: iextm_order = 1
101 integer :: last_tstep = -1
102 type(vector_list_t) :: interface_dof, interface_field
104 type(global_interpolation_settings_t) :: interpolation_settings
105 integer :: n_int_tot = 0
106 logical :: find_interface = .false.
107 logical :: setup = .false.
108 logical :: log = .false.
109
112 procedure(morph_overset_interface), nopass, pointer :: morph_interface => null()
113
114 contains
116 procedure, pass(this) :: init => overset_interface_vector_init
118 procedure, pass(this) :: init_from_components => &
121 procedure, pass(this) :: free => overset_interface_vector_free
123 procedure, pass(this) :: finalize => overset_interface_vector_finalize
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 => &
136 procedure, pass(this) :: update => overset_interface_update
137
139 procedure, pass(this), private :: build_masks_ => build_masks_
141 procedure, pass(this), private :: gather_interface_dofs_ => gather_interface_dofs_
143 procedure, pass(this), private :: setup_interpolator_ => setup_interpolator_
145 procedure, pass(this), private :: log_interface_error_ => log_interface_error_
146
148
149contains
150
154 subroutine overset_interface_vector_init(this, coef, json)
155 class(overset_interface_vector_t), intent(inout), target :: this
156 type(coef_t), target, intent(in) :: coef
157 type(json_file), intent(inout) ::json
158 real(kind=rp) :: tol, pad
159 logical :: log
160
162 call json_get_or_default(json, "interpolation.tolerance", &
163 tol, -1.0_rp)
164 call json_get_or_default(json, "interpolation.padding", &
165 pad, -1.0_rp)
166 call json_get_or_default(json, "order", this%iextm_order, 1)
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.")
169 end if
170 call json_get_or_default(json, "log", this%log, .false.)
171
172 call this%init_from_components(coef, tol, pad, log)
173
174 end subroutine overset_interface_vector_init
175
181 subroutine overset_interface_vector_init_from_components(this, coef, tol, pad, log)
182 class(overset_interface_vector_t), intent(inout), target :: this
183 type(coef_t), intent(in) :: coef
184 real(kind=rp), intent(in), optional :: tol, pad
185 logical, intent(in), optional :: log
186
188 call this%init_base(coef)
189
191 if (present(tol)) then
192 if (tol .gt. 0.0_rp) this%interpolation_settings%tolerance = tol
193 end if
194 if (present(pad)) then
195 if (pad .gt. 0.0_rp) this%interpolation_settings%padding = pad
196 end if
197 if (present(log)) then
198 this%log = log
199 end if
200
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")
204
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)
209
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')
214
217 if (neko_bcknd_device .eq. 1) then
218 ! copy
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())
222 ! synchronize
223 call this%x_dof%copy_from(device_to_host, sync = .false.)
224 call this%y_dof%copy_from(device_to_host, sync = .false.)
225 call this%z_dof%copy_from(device_to_host, sync = .true.)
226 else
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())
230 end if
231
233
237 class(overset_interface_vector_t), target, intent(inout) :: this
238
239 call this%bc_u%free()
240 call this%bc_v%free()
241 call this%bc_w%free()
242
243 call this%field_list%free()
244 call this%interface_dof%free()
245 call this%interface_field%free()
246
247 call this%x_dof%free()
248 call this%y_dof%free()
249 call this%z_dof%free()
250
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()
264
265 !if (associated(this%update_)) then
266 ! nullify(this%update_)
267 !end if
268 end subroutine overset_interface_vector_free
269
274 subroutine overset_interface_vector_apply_scalar(this, x, n, time, strong)
275 class(overset_interface_vector_t), intent(inout) :: this
276 integer, intent(in) :: n
277 real(kind=rp), intent(inout), dimension(n) :: x
278 type(time_state_t), intent(in), optional :: time
279 logical, intent(in), optional :: strong
280
281 call neko_error("overset_interface_vector cannot apply scalar BCs.&
282 & Use overset_interface_vector::apply_vector instead!")
283
285
289 subroutine overset_interface_vector_apply_scalar_dev(this, x_d, time, &
290 strong, strm)
291 class(overset_interface_vector_t), intent(inout), target :: this
292 type(c_ptr), intent(inout) :: x_d
293 type(time_state_t), intent(in), optional :: time
294 logical, intent(in), optional :: strong
295 type(c_ptr), intent(inout) :: strm
296
297 call neko_error("overset_interface_vector cannot apply scalar BCs.&
298 & Use overset_interface_vector::apply_vector instead!")
299
301
308 subroutine overset_interface_vector_apply_vector(this, x, y, z, n, time, &
309 strong)
310 class(overset_interface_vector_t), intent(inout) :: this
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
315 type(time_state_t), intent(in), optional :: time
316 logical, intent(in), optional :: strong
317 logical :: strong_
318
319 if (present(strong)) then
320 strong_ = strong
321 else
322 strong_ = .true.
323 end if
324
325 if (strong_) then
326
327 ! We can send any of the 3 bcs we have as argument, since they are all
328 ! the same boundary.
329 if (.not. this%updated) then
330 call this%update(time)
331 this%updated = .true.
332 end if
333
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))
338 end if
339
341
348 subroutine overset_interface_vector_apply_vector_dev(this, x_d, y_d, z_d, &
349 time, strong, strm)
350 class(overset_interface_vector_t), intent(inout), target :: this
351 type(c_ptr), intent(inout) :: x_d
352 type(c_ptr), intent(inout) :: y_d
353 type(c_ptr), intent(inout) :: z_d
354 type(time_state_t), intent(in), optional :: time
355 logical, intent(in), optional :: strong
356 type(c_ptr), intent(inout) :: strm
357 logical :: strong_
358
359 if (present(strong)) then
360 strong_ = strong
361 else
362 strong_ = .true.
363 end if
364
365 if (strong_) then
366 if (.not. this%updated) then
367 call this%update(time)
368 this%updated = .true.
369 end if
370
371 if (this%msk(0) .gt. 0) then
372 call device_masked_copy_0(x_d, this%bc_u%field_bc%x_d, &
373 this%bc_u%msk_d, this%bc_u%dof%size(), this%msk(0), strm) ! adperez: change the masks used here
374 call device_masked_copy_0(y_d, this%bc_v%field_bc%x_d, &
375 this%bc_v%msk_d, this%bc_v%dof%size(), this%msk(0), strm)
376 call device_masked_copy_0(z_d, this%bc_w%field_bc%x_d, &
377 this%bc_w%msk_d, this%bc_w%dof%size(), this%msk(0), strm)
378 end if
379 end if
380
382
384 subroutine overset_interface_vector_finalize(this, only_facets)
385 class(overset_interface_vector_t), target, intent(inout) :: this
386 logical, optional, intent(in) :: only_facets
387 logical :: only_facets_
388
389 if (present(only_facets)) then
390 only_facets_ = only_facets
391 else
392 only_facets_ = .false.
393 end if
394
396 call this%finalize_base(only_facets_)
397
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)
401
402 call this%bc_u%finalize(only_facets_)
403 call this%bc_v%finalize(only_facets_)
404 call this%bc_w%finalize(only_facets_)
405
407 call this%build_masks_()
408
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_()
414
416 call this%setup_interpolator_()
417
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')
422
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)
428
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)
433
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)
438
439 call mpi_allreduce(this%u_interface%size(), this%n_int_tot, 1, mpi_integer, &
440 mpi_sum, neko_global_comm)
441
442
444
446 subroutine overset_interface_update(this, time)
447 class(overset_interface_vector_t), intent(inout) :: this
448 type(time_state_t), intent(in) :: time
449 type(field_t), pointer :: u, v, w
450 type(iextm_time_scheme_t) :: time_scheme
451 integer :: nhist, ihist
452 real(kind=rp) :: iextm_coeffs(4)
453
454
456 call this%morph_interface(this%interface_dof, this%interface_field, &
457 this%interface_dof_mask, time, this%name, &
458 this%find_interface)
459
461 ! not implemented for now
462 ! if (substep .eq. 1) then
463 ! call this%extrapolate()
464 ! end if
465
467 if (this%find_interface) then
468
469 ! sync
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.)
473
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.
478
479 end if
480
482 u => neko_registry%get_field("u")
483 v => neko_registry%get_field("v")
484 w => neko_registry%get_field("w")
485
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.)
490
491 if (this%log) then
492 call this%log_interface_error_(u, v, w)
493 end if
494
495
497 if (time%tstep .ne. this%last_tstep) then
498 ! Update the last steps
499 this%last_tstep = time%tstep
500
501 ! Update the lag arrays with the value just got from the last step
502 call this%u_interface_lag%update()
503 call this%v_interface_lag%update()
504 call this%w_interface_lag%update()
505
506 ! Get the coefficients for the extrapolation
507 nhist = min(time%tstep, this%iextm_order)
508 call time_scheme%compute_coeffs(iextm_coeffs, time%dtlag, nhist)
509
510 ! Perfrom the extrapolation using the lag arrays
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))
514 do ihist = 2, nhist
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))
518 end do
519
520 end if
521
522
524 call vector_masked_scatter_copy(this%bc_u%field_bc%x(:,1,1,1), this%u_interface, &
525 this%interface_dof_mask, this%bc_u%dof%size())
526 call vector_masked_scatter_copy(this%bc_v%field_bc%x(:,1,1,1), this%v_interface, &
527 this%interface_dof_mask, this%bc_v%dof%size())
528 call vector_masked_scatter_copy(this%bc_w%field_bc%x(:,1,1,1), this%w_interface, &
529 this%interface_dof_mask, this%bc_w%dof%size())
530
531
532 end subroutine overset_interface_update
533
535 subroutine log_interface_error_(this, u, v, w)
536 class(overset_interface_vector_t), intent(inout) :: this
537 type(field_t), pointer, intent(in) :: u, v, w
538 real(kind=rp) :: u_int_norm, v_int_norm, w_int_norm
539 type(vector_t), pointer :: error
540 integer :: ind(1)
541 logical :: clear_scratch = .false.
542 character(len=256) :: log_buf
543
544 call neko_scratch_registry%request_vector(error, ind(1), this%u_interface%size(), &
545 clear_scratch)
546
548 call vector_masked_gather_copy(error, u%x(:,1,1,1), this%interface_dof_mask, &
549 this%dof%size())
550 call vector_add2s2(error, this%u_interface, -1.0_rp)
551 u_int_norm = sqrt(vector_glsc2(error, error)) / sqrt(real(this%n_int_tot, kind=rp))
552
553 call vector_masked_gather_copy(error, v%x(:,1,1,1), this%interface_dof_mask, &
554 this%dof%size())
555 call vector_add2s2(error, this%v_interface, -1.0_rp)
556 v_int_norm = sqrt(vector_glsc2(error, error)) / sqrt(real(this%n_int_tot, kind=rp))
557
558 call vector_masked_gather_copy(error, w%x(:,1,1,1), this%interface_dof_mask, &
559 this%dof%size())
560 call vector_add2s2(error, this%w_interface, -1.0_rp)
561 w_int_norm = sqrt(vector_glsc2(error, error)) / sqrt(real(this%n_int_tot, kind=rp))
562
563 call neko_scratch_registry%relinquish(ind)
564
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, ')'
569 call neko_log%message(log_buf)
570
571 end subroutine log_interface_error_
572
573 !===================
574 ! Helper subroutines
575 !===================
576
578 subroutine build_masks_(this)
579 class(overset_interface_vector_t), intent(inout) :: this
580 type(mask_t) :: temp_mask
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
585 type(stack_i4_t) :: stack
586
588 call this%interface_dof_mask%init(this%msk(1:this%msk(0)), this%msk(0))
589
592 lx = this%Xh%lx
593 ly = this%Xh%ly
594 lz = this%Xh%lz
595 allocate(found(this%msh%nelv))
596 found = .false.
597 !! Find sem elements that contain the boundary points
598 do i = 1, this%msk(0)
599 linear_idx = this%msk(i)
600 nonlinear_idx = nonlinear_index(linear_idx, lx, ly, lz)
601 found(nonlinear_idx(4)) = .true.
602 end do
603 !! fill the stack containing the gll indices
604 nelems = 0
605 call stack%init()
606 do e = 1, this%msh%nelv
607 if (found(e)) then
608 nelems = nelems + 1
609 do k = 1, this%Xh%lz
610 do j = 1, this%Xh%ly
611 do i = 1, this%Xh%lx
612 linear_idx = linear_index(i, j, k, e, lx, ly, lz)
613 call stack%push(linear_idx)
614 end do
615 end do
616 end do
617 end if
618 end do
619 deallocate(found)
620 call temp_mask%init(stack%array(), stack%size())
621 call stack%free()
622
624 call this%domain_element_mask%invert_mask(temp_mask, this%dof%size())
625
626 call temp_mask%free()
627
628 end subroutine build_masks_
629
631 subroutine gather_interface_dofs_(this)
632 class(overset_interface_vector_t), intent(inout) :: this
633
635 call vector_masked_gather_copy(this%x_interface_dof, &
636 this%dof%x(:,1,1,1), &
637 this%interface_dof_mask, &
638 this%dof%size())
639 call vector_masked_gather_copy(this%y_interface_dof, &
640 this%dof%y(:,1,1,1), &
641 this%interface_dof_mask, &
642 this%dof%size())
643 call vector_masked_gather_copy(this%z_interface_dof, &
644 this%dof%z(:,1,1,1), &
645 this%interface_dof_mask, &
646 this%dof%size())
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.)
651
652 end subroutine gather_interface_dofs_
653
654
656 subroutine setup_interpolator_(this)
657 class(overset_interface_vector_t), intent(inout) :: this
658
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)
665
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())
671
672 end subroutine setup_interpolator_
673
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
double real
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 list of bc_t.
Definition bc_list.f90:34
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 host_to_device
Definition device.F90:48
integer, parameter, public device_to_host
Definition device.F90:48
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
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 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.
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
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,...
Definition utils.f90:251
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_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.
Defines a vector.
Definition vector.f90:34
Base type for a boundary condition.
Definition bc.f90:62
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:49
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
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
Extension of the user defined dirichlet condition overset_interface
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...