Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
global_interpolation.f90
Go to the documentation of this file.
1! Copyright (c) 2020-2025, 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!
36 use num_types, only : rp, dp, xp
38 use space, only : space_t
39 use stack, only : stack_i4_t
40 use dofmap, only : dofmap_t
41 use logger, only : neko_log, log_size
43 use json_module, only : json_file
44 use utils, only : neko_error
54 use el_finder, only : el_finder_t
55 use pe_finder, only : pe_finder_t
56 use comm, only : neko_comm
57 use mpi_f08, only : mpi_sum, mpi_comm, mpi_comm_rank, &
58 mpi_comm_size, mpi_wtime, mpi_allreduce, mpi_in_place, mpi_integer, &
59 mpi_min, mpi_barrier, mpi_reduce_scatter_block, mpi_alltoall, &
60 mpi_isend, mpi_irecv
62 use vector, only : vector_t
64 use matrix, only : matrix_t
65 use math, only : copy, neko_eps
66 use mask, only : mask_t
67 use structs, only : array_ptr_t
68 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
69 implicit none
70 private
71
72 integer, public, parameter :: glob_map_size = 4096
73 real(kind=dp), public, parameter :: glob_interp_tol = neko_eps*1e3_dp
74 real(kind=dp), public, parameter :: glob_interp_pad = 1e-2_dp
75
81 real(kind=dp) :: tolerance = glob_interp_tol
83 real(kind=dp) :: padding = glob_interp_pad
85
87 type, public :: global_interpolation_t
89 type(vector_t) :: x
91 type(vector_t) :: y
93 type(vector_t) :: z
95 integer :: gdim
97 integer :: nelv
99 integer :: glb_nelv
101 type(mpi_comm) :: comm
103 integer :: pe_rank
105 integer :: pe_size
107 type(space_t) :: xh
110 integer :: n_points
112 real(kind=rp), allocatable :: xyz(:,:)
114 real(kind=rp), allocatable :: rst(:,:)
116 integer, allocatable :: pe_owner(:)
118 type(stack_i4_t), allocatable :: points_at_pe(:)
121 integer, allocatable :: el_owner0(:)
122 type(c_ptr) :: el_owner0_d = c_null_ptr
123
125 integer :: n_points_local
126 integer, allocatable :: el_owner0_local(:)
127 type(c_ptr) :: el_owner0_local_d = c_null_ptr
128 real(kind=rp), allocatable :: rst_local(:,:)
129 real(kind=rp), allocatable :: xyz_local(:,:)
131 type(local_interpolator_t) :: local_interp
134 logical :: all_points_local = .false.
136 real(kind=dp) :: tolerance = glob_interp_tol
138 real(kind=dp) :: padding = glob_interp_pad
139
143 integer, allocatable :: n_points_pe(:)
144 integer, allocatable :: n_points_offset_pe(:)
147 integer, allocatable :: n_points_pe_local(:)
148 integer, allocatable :: n_points_offset_pe_local(:)
151 class(pe_finder_t), allocatable :: pe_finder
153 class(el_finder_t), allocatable :: el_finder
155 type(legendre_rst_finder_t) :: rst_finder
160 type(vector_t) :: temp_local, temp
161 integer :: n_dof = -1
162 type(vector_t) :: masked_field
163
164 contains
167 procedure, pass(this) :: init_json_xyz => &
171 procedure, pass(this) :: init_json_dof => &
175 procedure, pass(this) :: init_xyz => global_interpolation_init_xyz
177 procedure, pass(this) :: init_dof => global_interpolation_init_dof
179 procedure, pass(this) :: free => global_interpolation_free
181 procedure, pass(this) :: free_points => global_interpolation_free_points
182 procedure, pass(this) :: free_points_local => &
184 procedure, pass(this) :: find_points_and_redist => &
189 procedure, pass(this) :: find_points_coords => &
191 procedure, pass(this) :: find_points_coords1d => &
194 procedure, pass(this) :: check_points => &
196 procedure, pass(this) :: find_points_xyz => global_interpolation_find_xyz
197 generic :: find_points => find_points_xyz, find_points_coords, &
198 find_points_coords1d
200 procedure, pass(this) :: evaluate => global_interpolation_evaluate
201 procedure, pass(this) :: evaluate_masked => &
203
205 generic :: init => init_dof, init_xyz, init_json_xyz, init_json_dof
206
208
209contains
210
224 subroutine global_interpolation_init_json_xyz(this, x, y, z, gdim, nelv, Xh, &
225 params_subdict, comm)
226 class(global_interpolation_t), target, intent(inout) :: this
227 real(kind=rp), intent(in) :: x(:)
228 real(kind=rp), intent(in) :: y(:)
229 real(kind=rp), intent(in) :: z(:)
230 integer, intent(in) :: gdim
231 integer, intent(in) :: nelv
232 type(space_t), intent(in) :: Xh
233 type(json_file), intent(inout) :: params_subdict
234 type(mpi_comm), intent(in), optional :: comm
235
236 real(kind=dp) :: tol, pad
237
238 call json_get_or_lookup_or_default(params_subdict, 'tolerance', &
239 tol, glob_interp_tol)
240 call json_get_or_lookup_or_default(params_subdict, 'padding', &
241 pad, glob_interp_pad)
242
243 call this%init_xyz(x, y, z, gdim, nelv, xh, comm = comm, tol = tol, &
244 pad = pad)
245
247
255 subroutine global_interpolation_init_json_dof(this, dof, params_subdict, &
256 comm, mask)
257 class(global_interpolation_t), target, intent(inout) :: this
258 type(dofmap_t) :: dof
259 type(json_file), intent(inout) :: params_subdict
260 type(mpi_comm), optional, intent(in) :: comm
261 type(mask_t), intent(in), optional :: mask
262
263 real(kind=dp) :: tol, pad
264
265 call json_get_or_lookup_or_default(params_subdict, 'tolerance', &
266 tol, glob_interp_tol)
267 call json_get_or_lookup_or_default(params_subdict, 'padding', &
268 pad, glob_interp_pad)
269
270 call this%init_dof(dof, comm = comm, tol = tol, pad = pad, mask = mask)
271
273
280 subroutine global_interpolation_init_dof(this, dof, comm, tol, pad, mask)
281 class(global_interpolation_t), target, intent(inout) :: this
282 type(dofmap_t) :: dof
283 type(mpi_comm), optional, intent(in) :: comm
284 real(kind=dp), optional :: tol
285 real(kind=dp), optional :: pad
286 type(mask_t), intent(in), optional :: mask
287
288 integer :: temp_nelv
289
290 ! Store the number of dofs
291 this%n_dof = dof%size()
292 ! NOTE: Passing dof%x(:,1,1,1), etc in init_xyz passes down the entire
293 ! dof%x array and not a slice. It is done this way for
294 ! to get the right dimension (see global_interpolation_init_xyz).
295 if (.not. present(mask)) then
296 call this%init_xyz(dof%x(:,1,1,1), dof%y(:,1,1,1), dof%z(:,1,1,1), &
297 dof%msh%gdim, dof%msh%nelv, dof%Xh, comm = comm, &
298 tol = tol, pad = pad)
299 else
300
301 ! Initialize a helper field with the size of the mask
302 call this%masked_field%init(mask%size())
303 ! Verify that the mask size is compatible with the dofmap
304 temp_nelv = mask%size() / (dof%Xh%lx*dof%Xh%ly*dof%Xh%lz)
305 if (mod(mask%size(), dof%Xh%lx*dof%Xh%ly*dof%Xh%lz) /= 0) then
306 call neko_error("Mask size must be a multiple of the number of" // &
307 " elements in the mesh.")
308 end if
309 ! Initialize with the masked coordinates
310 call this%init_xyz(dof%x(mask%get(),1,1,1), dof%y(mask%get(),1,1,1), &
311 dof%z(mask%get(),1,1,1), dof%msh%gdim, temp_nelv, dof%Xh, &
312 comm = comm, tol = tol, pad = pad)
313 end if
314
315 end subroutine global_interpolation_init_dof
316
327 subroutine global_interpolation_init_xyz(this, x, y, z, gdim, nelv, Xh, &
328 comm, tol, pad)
329 class(global_interpolation_t), target, intent(inout) :: this
330 real(kind=rp), intent(in) :: x(:)
331 real(kind=rp), intent(in) :: y(:)
332 real(kind=rp), intent(in) :: z(:)
333 integer, intent(in) :: gdim
334 integer, intent(in) :: nelv
335 type(space_t), intent(in) :: Xh
336 type(mpi_comm), intent(in), optional :: comm
337 real(kind=dp), intent(in), optional :: tol
338 real(kind=dp), intent(in), optional :: pad
339
340 integer :: lx, ly, lz, ierr, i, n
341 character(len=8000) :: log_buf
342 !real(kind=dp) :: padding ! <-- note that this padding is in dp
343 real(kind=rp) :: time1, time_start
344 character(len=255) :: mode_str
345 integer :: boxdim, envvar_len
346
347 call neko_log%section('Global Interpolation')
348 call neko_log%message('Initializing global interpolation')
349
350 call this%free()
351
352 ! Set communicator
353 if (present(comm)) then
354 this%comm = comm
355 else
356 this%comm = neko_comm
357 end if
358
359 ! Set point search parameters
360 this%padding = glob_interp_pad
361 if (present(pad)) this%padding = pad
362
363 this%tolerance = glob_interp_tol
364 if (present(tol)) this%tolerance = tol
365
366 write(log_buf, '(A,E15.7)') &
367 'Tolerance: ', this%tolerance
368 call neko_log%message(log_buf)
369 write(log_buf, '(A,E15.7)') &
370 'Padding : ', this%padding
371 call neko_log%message(log_buf)
372
373 time_start = mpi_wtime()
374 call mpi_barrier(this%comm)
375
376 call mpi_comm_rank(this%comm, this%pe_rank, ierr)
377 call mpi_comm_size(this%comm, this%pe_size, ierr)
378
379 this%gdim = gdim
380 this%nelv = nelv
381
382 call mpi_allreduce(nelv, this%glb_nelv, 1, mpi_integer, &
383 mpi_sum, this%comm, ierr)
384 lx = xh%lx
385 ly = xh%ly
386 lz = xh%lz
387 n = nelv * lx*ly*lz
388 call this%x%init(n)
389 call this%y%init(n)
390 call this%z%init(n)
391 call copy(this%x%x, x, n)
392 call this%x%copy_from(host_to_device,.false.)
393 call copy(this%y%x, y, n)
394 call this%y%copy_from(host_to_device,.false.)
395 call copy(this%z%x, z, n)
396 call this%z%copy_from(host_to_device,.false.)
397 call this%Xh%init(xh%t, lx, ly, lz)
398
399 ! Initialize n_dof if it was not started with a dof-based constructor
400 if (this%n_dof == -1) then
401 this%n_dof = n
402 end if
403
404
405 call get_environment_variable("NEKO_GLOBAL_INTERP_EL_FINDER", &
406 mode_str, envvar_len)
407
408 if (envvar_len .gt. 0) then
409 if (mode_str(1:envvar_len) == 'AABB') then
410 allocate(aabb_el_finder_t :: this%el_finder)
411 end if
412 end if
413 call get_environment_variable("NEKO_GLOBAL_INTERP_PE_FINDER", &
414 mode_str, envvar_len)
415
416 if (envvar_len .gt. 0) then
417 if (mode_str(1:envvar_len) == 'AABB') then
418 allocate(aabb_pe_finder_t :: this%pe_finder)
419 end if
420 end if
421
422 if (.not. allocated(this%el_finder)) then
423 allocate(cartesian_el_finder_t :: this%el_finder)
424 end if
425 if (.not. allocated(this%pe_finder)) then
426 allocate(cartesian_pe_finder_t :: this%pe_finder)
427 end if
428 select type (el_find => this%el_finder)
429 type is (aabb_el_finder_t)
430 call neko_log%message('Using AABB element finder')
431 call el_find%init(x, y, z, nelv, xh, this%padding)
432 type is (cartesian_el_finder_t)
433 call neko_log%message('Using Cartesian element finder')
434 boxdim = max(lx*int(real(nelv, xp)**(1.0_xp / 3.0_xp)), 2)
435 boxdim = min(boxdim, 300)
436 call el_find%init(x, y, z, nelv, xh, boxdim, this%padding)
437 class default
438 call neko_error('Unknown element finder type')
439 end select
440
441 select type (pe_find => this%pe_finder)
442 type is (aabb_pe_finder_t)
443 call neko_log%message('Using AABB PE finder')
444 call pe_find%init(this%x%x, this%y%x, this%z%x, &
445 nelv, xh, this%comm, this%padding)
446 type is (cartesian_pe_finder_t)
447 call neko_log%message('Using Cartesian PE finder')
448 boxdim = lx*int(real(this%glb_nelv, xp)**(1.0_xp / 3.0_xp))
449 boxdim = max(boxdim, 32)
450 boxdim = min(boxdim, &
451 int(8.0_xp*(30000.0_xp * this%pe_size)**(1.0_xp / 3.0_xp)))
452 call pe_find%init(this%x%x, this%y%x, this%z%x, &
453 nelv, xh, this%comm, boxdim, this%padding)
454 class default
455 call neko_error('Unknown PE finder type')
456 end select
457
458 call this%rst_finder%init(this%x%x, this%y%x, this%z%x, nelv, xh, &
459 this%tolerance)
460 if (allocated(this%n_points_pe)) deallocate(this%n_points_pe)
461 if (allocated(this%n_points_pe_local)) deallocate(this%n_points_pe_local)
462 if (allocated(this%n_points_offset_pe_local)) &
463 deallocate(this%n_points_offset_pe_local)
464 if (allocated(this%n_points_offset_pe)) deallocate(this%n_points_offset_pe)
465 allocate(this%n_points_pe(0:(this%pe_size-1)))
466 allocate(this%n_points_offset_pe(0:(this%pe_size-1)))
467 allocate(this%n_points_pe_local(0:(this%pe_size-1)))
468 allocate(this%n_points_offset_pe_local(0:(this%pe_size-1)))
469 allocate(this%points_at_pe(0:(this%pe_size-1)))
470 do i = 0, this%pe_size-1
471 call this%points_at_pe(i)%init()
472 end do
473 call mpi_barrier(this%comm)
474 time1 = mpi_wtime()
475 write(log_buf, '(A,E15.7)') &
476 'Global interpolation initialized (s):', time1-time_start
477 call neko_log%message(log_buf)
478 call neko_log%end_section()
479 end subroutine global_interpolation_init_xyz
480
481
484 class(global_interpolation_t), target, intent(inout) :: this
485 integer :: i
486
487 call this%x%free()
488 call this%y%free()
489 call this%z%free()
490 call this%Xh%free()
491
492 this%nelv = 0
493 this%gdim = 0
494
495 call this%free_points()
496 call this%free_points_local()
497 call this%local_interp%free()
498 if (allocated(this%el_finder)) then
499 call this%el_finder%free()
500 deallocate(this%el_finder)
501 end if
502 if (allocated(this%pe_finder)) then
503 call this%pe_finder%free()
504 deallocate(this%pe_finder)
505 end if
506 call this%rst_finder%free()
507
508 call this%temp_local%free()
509 call this%temp%free()
510 if (allocated(this%points_at_pe)) then
511 do i = 0, this%pe_size-1
512 call this%points_at_pe(i)%free()
513 end do
514 deallocate(this%points_at_pe)
515 end if
516 if (allocated(this%n_points_pe)) deallocate(this%n_points_pe)
517 if (allocated(this%n_points_pe_local)) deallocate(this%n_points_pe_local)
518 if (allocated(this%n_points_offset_pe_local)) &
519 deallocate(this%n_points_offset_pe_local)
520 if (allocated(this%n_points_offset_pe)) deallocate(this%n_points_offset_pe)
521
522
523
524 end subroutine global_interpolation_free
525
528 class(global_interpolation_t), target, intent(inout) :: this
529
530 this%n_points = 0
531 this%all_points_local = .false.
532
533 if (allocated(this%xyz)) deallocate(this%xyz)
534 if (allocated(this%rst)) deallocate(this%rst)
535 if (allocated(this%pe_owner)) deallocate(this%pe_owner)
536 if (allocated(this%el_owner0)) then
537 if (neko_bcknd_device .eq. 1) then
538 call device_unmap(this%el_owner0, this%el_owner0_d)
539 end if
540 deallocate(this%el_owner0)
541 end if
542
543 call this%glb_intrp_comm%free()
544
545
547
548
550 class(global_interpolation_t), target, intent(inout) :: this
551
552 this%n_points_local = 0
553 this%all_points_local = .false.
554
555 if (allocated(this%xyz_local)) deallocate(this%xyz_local)
556 if (allocated(this%rst_local)) deallocate(this%rst_local)
557 if (allocated(this%el_owner0_local)) then
558 if (neko_bcknd_device .eq. 1) then
559 call device_unmap(this%el_owner0_local, this%el_owner0_local_d)
560 end if
561 deallocate(this%el_owner0_local)
562 end if
563
565
566
569 class(global_interpolation_t), target, intent(inout) :: this
570 character(len=8000) :: log_buf
571 type(vector_t) :: x_t
572 type(vector_t) :: y_t
573 type(vector_t) :: z_t
574 type(matrix_t) :: rst_local_cand
575 type(vector_t) :: resx
576 type(vector_t) :: resy
577 type(vector_t) :: resz
578 type(c_ptr) :: el_cands_d
579 type(matrix_t) :: res
580 integer :: i, j, stupid_intent
581 type(stack_i4_t), target :: all_el_candidates
582 integer, allocatable :: n_el_cands(:)
583 integer, contiguous, pointer :: el_cands(:), point_ids(:), send_recv(:)
584 real(kind=rp), allocatable :: res_results(:,:)
585 real(kind=rp), allocatable :: rst_results(:,:)
586 integer, allocatable :: el_owner_results(:)
587 integer :: ierr, ii, n_point_cand, n_glb_point_cand, point_id, rank
588 real(kind=rp) :: time1, time2, time_start
589 !Temp stuff for glb_intrp_comm
590 type(stack_i4_t) :: send_pe, recv_pe
591 type(glb_intrp_comm_t) :: glb_intrp_find, glb_intrp_find_back
592 type(stack_i4_t) :: send_pe_find, recv_pe_find
593
594 el_cands_d = c_null_ptr
595 call neko_log%section('Global Interpolation')
596 call glb_intrp_find%init_dofs(this%pe_size)
597 call send_pe_find%init()
598 call recv_pe_find%init()
599 call mpi_barrier(this%comm)
600 time_start = mpi_wtime()
601 write(log_buf, '(A)') 'Global interpolation, finding points'
602 call neko_log%message(log_buf)
603 ! Find pe candidates that the points i want may be at
604 ! Add number to n_points_pe_local
605 !Working arrays
606 this%n_points_pe = 0
607 call this%pe_finder%find_batch(this%xyz, this%n_points, &
608 this%points_at_pe, this%n_points_pe)
609 call mpi_barrier(this%comm)
610 time1 = mpi_wtime()
611 write(log_buf, '(A,E15.7)') &
612 'Found PE candidates time since start of findpts (s):', &
613 time1 - time_start
614 call neko_log%message(log_buf)
615
616 !Send number of points I want to candidates
617 ! n_points_local -> how many points might be at this rank
618 ! n_points_pe_local -> how many points local on this rank that other pes
619 ! might want
620 this%n_points_pe_local = 0
621 this%n_points_local = 0
622 call mpi_reduce_scatter_block(this%n_points_pe, this%n_points_local, &
623 1, mpi_integer, mpi_sum, this%comm, ierr)
624 call mpi_alltoall(this%n_points_pe, 1, mpi_integer,&
625 this%n_points_pe_local, 1, mpi_integer, this%comm, ierr)
626
627 !Set up offset arrays
628 this%n_points_offset_pe_local(0) = 0
629 this%n_points_offset_pe(0) = 0
630 do i = 1, (this%pe_size - 1)
631 this%n_points_offset_pe_local(i) = this%n_points_pe_local(i-1)&
632 + this%n_points_offset_pe_local(i-1)
633 this%n_points_offset_pe(i) = this%n_points_pe(i-1)&
634 + this%n_points_offset_pe(i-1)
635 end do
636 do i = 0, (this%pe_size-1)
637 if (this%n_points_pe(i) .gt. 0) then
638 call send_pe_find%push(i)
639 point_ids => this%points_at_pe(i)%array()
640 do j = 1, this%n_points_pe(i)
641 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+1)
642 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+2)
643 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+3)
644 end do
645 end if
646 if (this%n_points_pe_local(i) .gt. 0) then
647 call recv_pe_find%push(i)
648 do j = 1, this%n_points_pe_local(i)
649 call glb_intrp_find%recv_dof(i)%push(3*(j + &
650 this%n_points_offset_pe_local(i) - 1) + 1)
651 call glb_intrp_find%recv_dof(i)%push(3*(j + &
652 this%n_points_offset_pe_local(i) - 1) + 2)
653 call glb_intrp_find%recv_dof(i)%push(3*(j + &
654 this%n_points_offset_pe_local(i) - 1) + 3)
655 end do
656 end if
657 end do
658
659
660
661 call glb_intrp_find%init(send_pe_find, recv_pe_find, this%comm)
662
663 call glb_intrp_find_back%init_dofs(this%pe_size)
664 ii = 0
665 do i = 0, (this%pe_size-1)
666 send_recv => glb_intrp_find%recv_dof(i)%array()
667 do j = 1, glb_intrp_find%recv_dof(i)%size()
668 call glb_intrp_find_back%send_dof(i)%push(send_recv(j))
669 end do
670 send_recv => glb_intrp_find%send_dof(i)%array()
671 do j = 1, glb_intrp_find%send_dof(i)%size()
672 ii = ii + 1
673 call glb_intrp_find_back%recv_dof(i)%push(ii)
674 end do
675 end do
676
677 call glb_intrp_find_back%init(recv_pe_find, send_pe_find, this%comm)
678
679
680 if (allocated(this%xyz_local)) then
681 deallocate(this%xyz_local)
682 end if
683 allocate(this%xyz_local(3, this%n_points_local))
684 call glb_intrp_find%sendrecv(this%xyz, this%xyz_local, this%n_points*3, &
685 this%n_points_local*3)
686
687 call mpi_barrier(this%comm)
688 time1 = mpi_wtime()
689 write(log_buf, '(A,E15.7)') &
690 'Sent to points to PE candidates, time since start of ' &
691 // 'find_points (s):', time1 - time_start
692 call neko_log%message(log_buf)
693
694 !Okay, now we need to find the rst...
695 call all_el_candidates%init()
696
697 if (allocated(n_el_cands)) then
698 deallocate(n_el_cands)
699 end if
700
701 allocate(n_el_cands(this%n_points_local))
703 call this%el_finder%find_batch(this%xyz_local, this%n_points_local, &
704 all_el_candidates, n_el_cands)
705
706 n_point_cand = all_el_candidates%size()
707 if (n_point_cand .gt. 1e8) then
708 print *, 'Warning, many point candidates on rank', this%pe_rank, &
709 'cands:', n_point_cand, &
710 'Consider increasing number of ranks'
711 end if
712 call x_t%init(n_point_cand)
713 call y_t%init(n_point_cand)
714 call z_t%init(n_point_cand)
715 ii = 0
717 do i = 1 , this%n_points_local
718 do j = 1, n_el_cands(i)
719 ii = ii + 1
720 x_t%x(ii) = this%xyz_local(1,i)
721 y_t%x(ii) = this%xyz_local(2,i)
722 z_t%x(ii) = this%xyz_local(3,i)
723 end do
724 end do
725
726 call mpi_barrier(this%comm)
727 time1 = mpi_wtime()
728 write(log_buf, '(A,E15.7)') &
729 'Element candidates found, now time for finding rst, time ' // &
730 'since start of find_points (s):', time1 - time_start
731 call neko_log%message(log_buf)
732 call rst_local_cand%init(3, n_point_cand)
733 call resx%init(n_point_cand)
734 call resy%init(n_point_cand)
735 call resz%init(n_point_cand)
736
737 ! Find rst within all element candidates for target xyz (x_t, y_t, z_t)
738 call mpi_barrier(this%comm)
739 time1 = mpi_wtime()
740 el_cands => all_el_candidates%array()
741 if (neko_bcknd_device .eq. 1) then
742 call x_t%copy_from(host_to_device,.false.)
743
744 call y_t%copy_from(host_to_device,.false.)
745
746 call z_t%copy_from(host_to_device,.false.)
747 call device_map(el_cands, el_cands_d, n_point_cand)
748 call device_memcpy(el_cands, el_cands_d, n_point_cand, &
749 host_to_device, .true.)
750 end if
751
752 call this%rst_finder%find(rst_local_cand, &
753 x_t, y_t, z_t, &
754 el_cands, n_point_cand, &
755 resx, resy, resz)
756 if (neko_bcknd_device .eq. 1) then
757 call rst_local_cand%copy_from(device_to_host,.false.)
758 call resx%copy_from(device_to_host,.false.)
759 call resy%copy_from(device_to_host,.false.)
760 call resz%copy_from(device_to_host,.true.)
761 call device_unmap(el_cands, el_cands_d)
762 end if
763 call mpi_barrier(this%comm)
764
765 time2 = mpi_wtime()
766 write(log_buf, '(A,E15.7)') &
767 'Found rst with Newton iteration, time (s):', time2-time1
768 call neko_log%message(log_buf)
769
770 write(log_buf, '(A)') &
771 'Checking validity of points and choosing best candidates.'
772 call neko_log%message(log_buf)
773 call mpi_barrier(this%comm, ierr)
774
775 if (allocated(this%rst_local)) deallocate(this%rst_local)
776 if (allocated(this%el_owner0_local)) deallocate(this%el_owner0_local)
777 allocate(this%rst_local(3, this%n_points_local))
778 allocate(this%el_owner0_local(this%n_points_local))
779 ! Choose the best candidate at this rank
780 ii = 0
781 do i = 1 , this%n_points_local
782 this%xyz_local(1,i) = 10.0
783 this%xyz_local(2,i) = 10.0
784 this%xyz_local(3,i) = 10.0
785 this%rst_local(1,i) = 10.0
786 this%rst_local(2,i) = 10.0
787 this%rst_local(3,i) = 10.0
788 this%el_owner0_local(i) = -1
789 do j = 1, n_el_cands(i)
790 ii = ii + 1
791 if (rst_cmp(this%rst_local(:, i), rst_local_cand%x(:, ii), &
792 this%xyz_local(:, i), [resx%x(ii), resy%x(ii), resz%x(ii)], &
793 this%padding)) then
794 this%rst_local(1, i) = rst_local_cand%x(1, ii)
795 this%rst_local(2, i) = rst_local_cand%x(2, ii)
796
797 this%rst_local(3, i) = rst_local_cand%x(3, ii)
798 this%xyz_local(1,i) = resx%x(ii)
799 this%xyz_local(2,i) = resy%x(ii)
800 this%xyz_local(3,i) = resz%x(ii)
801 this%el_owner0_local(i) = el_cands(ii)
802 end if
803 ! if (this%pe_rank .eq. 0) print *,i, this%rst_local(:,i), &
804 ! this%xyz_local(:,i), this%el_owner0_local(i)
805 end do
806 end do
807 call res%init(3, this%n_points)
808 n_glb_point_cand = sum(this%n_points_pe)
809 if (allocated(rst_results)) deallocate(rst_results)
810 if (allocated(res_results)) deallocate(res_results)
811 if (allocated(el_owner_results)) deallocate(el_owner_results)
812 allocate(rst_results(3, n_glb_point_cand))
813 allocate(res_results(3, n_glb_point_cand))
814 allocate(el_owner_results(n_glb_point_cand))
815 res = 1e2_rp
816 this%rst = 1e2
817 this%pe_owner = -1
818 this%el_owner0 = -1
819 call glb_intrp_find_back%sendrecv(this%xyz_local, res_results, &
820 this%n_points_local*3, n_glb_point_cand*3)
821 call glb_intrp_find_back%sendrecv(this%rst_local, rst_results, &
822 this%n_points_local*3, n_glb_point_cand*3)
823 do i = 1, size(glb_intrp_find_back%send_pe)
824 rank = glb_intrp_find_back%send_pe(i)
825 call mpi_isend(this%el_owner0_local( &
826 this%n_points_offset_pe_local(rank) + 1), &
827 this%n_points_pe_local(rank), &
828 mpi_integer, rank, 0, &
829 this%comm, glb_intrp_find_back%send_buf(i)%request, ierr)
830 glb_intrp_find_back%send_buf(i)%flag = .false.
831 end do
832 do i = 1, size(glb_intrp_find_back%recv_pe)
833 rank = glb_intrp_find_back%recv_pe(i)
834 call mpi_irecv(el_owner_results(this%n_points_offset_pe(rank)+1),&
835 this%n_points_pe(rank), &
836 mpi_integer, rank, 0, &
837 this%comm, glb_intrp_find_back%recv_buf(i)%request, ierr)
838 glb_intrp_find_back%recv_buf(i)%flag = .false.
839 end do
840 call glb_intrp_find_back%nbwait_no_op()
841 ii = 0
842 do i = 1, size(glb_intrp_find_back%recv_pe)
843 point_ids => this%points_at_pe(glb_intrp_find_back%recv_pe(i))%array()
844 do j = 1, this%n_points_pe(glb_intrp_find_back%recv_pe(i))
845 point_id = point_ids(j)
846 ii = ii + 1
847 if (rst_cmp(this%rst(:, point_id), rst_results(:, ii), &
848 res%x(:, point_id), res_results(:, ii), this%padding) .or. &
849 this%pe_owner(point_ids(j)) .eq. -1 ) then
850 this%rst(:, point_ids(j)) = rst_results(:, ii)
851 res%x(:, point_ids(j)) = res_results(:, ii)
852 this%pe_owner(point_ids(j)) = glb_intrp_find_back%recv_pe(i)
853 this%el_owner0(point_ids(j)) = el_owner_results(ii)
854 end if
855 ! if (this%pe_rank .eq. 0) print *,point_id, &
856 !this%rst(:,point_ids(j)),res%x(:,point_ids(j)),
857 ! this%el_owner0(point_ids(j))
858 end do
859 end do
860
861 !OK, now I know the correct rst values of the points I want We now
862 !send the correct rsts to the correct rank (so a point only
863 !belongs to one rank)
864 do i = 0, this%pe_size-1
865 call this%points_at_pe(i)%clear()
866 this%n_points_pe(i) = 0
867 end do
868
869 do i = 1, this%n_points
870 stupid_intent = i
871 if (this%pe_owner(i) .eq. -1 .or. this%el_owner0(i) .eq. -1) then
872 print *, 'No owning rank found for',&
873 ' point ', stupid_intent, ' with coords', this%xyz(:,i), &
874 ' Interpolation will always yield 0.0. Try increase padding.'
875 else
876 call this%points_at_pe(this%pe_owner(i))%push(stupid_intent)
877 this%n_points_pe(this%pe_owner(i)) = &
878 this%n_points_pe(this%pe_owner(i)) + 1
879 end if
880 end do
881 call mpi_reduce_scatter_block(this%n_points_pe, this%n_points_local, 1, &
882 mpi_integer, mpi_sum, this%comm, ierr)
883 call mpi_alltoall(this%n_points_pe, 1, mpi_integer, &
884 this%n_points_pe_local, 1, mpi_integer, this%comm, ierr)
885 this%n_points_offset_pe_local(0) = 0
886 this%n_points_offset_pe(0) = 0
887 do i = 1, (this%pe_size - 1)
888 this%n_points_offset_pe_local(i) = this%n_points_pe_local(i-1)&
889 + this%n_points_offset_pe_local(i-1)
890 this%n_points_offset_pe(i) = this%n_points_pe(i-1)&
891 + this%n_points_offset_pe(i-1)
892 end do
893 call send_pe_find%free()
894 call recv_pe_find%free()
895 call glb_intrp_find%free()
896 call send_pe_find%init()
897 call recv_pe_find%init()
898 call glb_intrp_find%init_dofs(this%pe_size)
899 !setup comm to send xyz and rst to chosen ranks
900 do i = 0, (this%pe_size-1)
901 if (this%n_points_pe(i) .gt. 0) then
902 call send_pe_find%push(i)
903 point_ids => this%points_at_pe(i)%array()
904 do j = 1, this%n_points_pe(i)
905 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 1)
906 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 2)
907 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 3)
908 end do
909 end if
910 if (this%n_points_pe_local(i) .gt. 0) then
911 call recv_pe_find%push(i)
912 do j = 1, this%n_points_pe_local(i)
913 call glb_intrp_find%recv_dof(i)%push(3*(j + &
914 this%n_points_offset_pe_local(i) - 1) + 1)
915 call glb_intrp_find%recv_dof(i)%push(3*(j + &
916 this%n_points_offset_pe_local(i) - 1) + 2)
917 call glb_intrp_find%recv_dof(i)%push(3*(j + &
918 this%n_points_offset_pe_local(i) - 1) + 3)
919 end do
920 end if
921 end do
922
923
924 call glb_intrp_find%init(send_pe_find, recv_pe_find, this%comm)
925 call glb_intrp_find%sendrecv(this%xyz, this%xyz_local, this%n_points*3, &
926 this%n_points_local*3)
927 call glb_intrp_find%sendrecv(this%rst, this%rst_local, this%n_points*3, &
928 this%n_points_local*3)
929 ii = 0
930 do i = 1, size(glb_intrp_find%send_pe)
931 rank = glb_intrp_find%send_pe(i)
932 point_ids => this%points_at_pe(rank)%array()
933 do j = 1, this%n_points_pe(rank)
934 ii = ii + 1
935 el_owner_results(ii) = this%el_owner0(point_ids(j))
936 end do
937 call mpi_isend(el_owner_results(this%n_points_offset_pe(rank) + 1),&
938 this%n_points_pe(rank), &
939 mpi_integer, rank, 0, &
940 this%comm, glb_intrp_find%send_buf(i)%request, ierr)
941 glb_intrp_find%send_buf(i)%flag = .false.
942 end do
943 do i = 1, size(glb_intrp_find%recv_pe)
944 rank = glb_intrp_find%recv_pe(i)
945 call mpi_irecv(this%el_owner0_local( &
946 this%n_points_offset_pe_local(rank) + 1), &
947 this%n_points_pe_local(rank), &
948 mpi_integer, rank, 0, &
949 this%comm, glb_intrp_find%recv_buf(i)%request, ierr)
950 glb_intrp_find%recv_buf(i)%flag = .false.
951 end do
952 call glb_intrp_find%nbwait_no_op()
953
954 call glb_intrp_find%free()
955
956 !Set up final way of doing communication
957 call send_pe%init()
958 call recv_pe%init()
959 call this%glb_intrp_comm%init_dofs(this%pe_size)
960 do i = 0, (this%pe_size-1)
961 if (this%n_points_pe(i) .gt. 0) then
962 call recv_pe%push(i)
963 point_ids => this%points_at_pe(i)%array()
964 do j = 1, this%n_points_pe(i)
965 call this%glb_intrp_comm%recv_dof(i)%push(point_ids(j))
966 end do
967 end if
968 if (this%n_points_pe_local(i) .gt. 0) then
969 call send_pe%push(i)
970 do j = 1, this%n_points_pe_local(i)
971 call this%glb_intrp_comm%send_dof(i)%push(j + &
972 this%n_points_offset_pe_local(i))
973 end do
974 end if
975 end do
976 call this%glb_intrp_comm%init(send_pe, recv_pe, this%comm)
977
978 !Initialize working arrays for evaluation
979 call this%temp_local%init(this%n_points_local)
980 call this%temp%init(this%n_points)
981
982 !Initialize interpolator for local interpolation
983 call this%local_interp%init(this%Xh, this%rst_local, &
984 this%n_points_local)
985
986
987 if (neko_bcknd_device .eq. 1) then
988 call device_memcpy(this%el_owner0, this%el_owner0_d, &
989 this%n_points, host_to_device, sync = .true.)
990 call device_map(this%el_owner0_local, this%el_owner0_local_d, &
991 this%n_points_local)
992 call device_memcpy(this%el_owner0_local, this%el_owner0_local_d, &
993 this%n_points_local, host_to_device, sync = .true.)
994 end if
995
996 call this%check_points(this%x%x, this%y%x, this%z%x)
997
998 !Free stuff
999 call send_pe%free()
1000 call recv_pe%free()
1001 call glb_intrp_find_back%free()
1002 call send_pe_find%free()
1003 call recv_pe_find%free()
1004 call x_t%free()
1005 call y_t%free()
1006 call z_t%free()
1007 call rst_local_cand%free()
1008 call resx%free()
1009 call resy%free()
1010 call resz%free()
1011 call res%free()
1012 call all_el_candidates%free()
1013
1014 if (allocated(n_el_cands)) deallocate(n_el_cands)
1015 if (allocated(rst_results)) deallocate(rst_results)
1016 if (allocated(res_results)) deallocate(res_results)
1017 if (allocated(el_owner_results)) deallocate(el_owner_results)
1018 call mpi_barrier(this%comm, ierr)
1019 time2 = mpi_wtime()
1020 write(log_buf, '(A,E15.7)') 'Global interpolation find points ' // &
1021 'done, time (s):', time2-time_start
1022 call neko_log%message(log_buf)
1023 call neko_log%end_section()
1024 call neko_log%newline()
1026
1030 subroutine global_interpolation_check_points(this, x, y, z)
1031 class(global_interpolation_t), target, intent(inout) :: this
1032 real(kind=rp), intent(inout) :: x(:)
1033 real(kind=rp), intent(inout) :: y(:)
1034 real(kind=rp), intent(inout) :: z(:)
1035 integer :: i, j
1036 character(len=8000) :: log_buf
1037 real(kind=rp) :: xdiff, ydiff, zdiff
1038 logical :: isdiff
1039 type(vector_t) :: x_check, y_check, z_check
1040
1041 call x_check%init(this%n_points)
1042 call y_check%init(this%n_points)
1043 call z_check%init(this%n_points)
1044 call this%evaluate(x_check%x, x, on_host = .true.)
1045 call this%evaluate(y_check%x, y, on_host = .true.)
1046 call this%evaluate(z_check%x, z, on_host = .true.)
1047 write(log_buf, '(A)') 'Checking validity of points.'
1048 call neko_log%message(log_buf)
1049 j = 0
1050 do i = 1 , this%n_points
1051 ! Check validity of points
1052 isdiff = .false.
1053 xdiff = x_check%x(i)-this%xyz(1,i)
1054 ydiff = y_check%x(i)-this%xyz(2,i)
1055 zdiff = z_check%x(i)-this%xyz(3,i)
1056 isdiff = norm2(real([xdiff, ydiff, zdiff], xp)) > this%tolerance
1057 if (isdiff) then
1058 write(*, *) 'Point ', i, 'at rank ', this%pe_rank, &
1059 'with coordinates: ', &
1060 this%xyz(1, i), this%xyz(2, i), this%xyz(3, i), &
1061 'Differ from interpolated coords: ', &
1062 x_check%x(i), y_check%x(i), z_check%x(i), &
1063 'Actual difference: ', &
1064 xdiff, ydiff, zdiff, norm2(real([xdiff, ydiff, zdiff], xp)), &
1065 'Process, element: ', &
1066 this%pe_owner(i), this%el_owner0(i)+1, &
1067 'Calculated rst: ', &
1068 this%rst(1,i), this%rst(2,i), this%rst(3,i)
1069 j = j + 1
1070 end if
1071 end do
1072 call x_check%free()
1073 call y_check%free()
1074 call z_check%free()
1075
1077
1085 subroutine global_interpolation_find_coords(this, x, y, z, n_points)
1086 class(global_interpolation_t), intent(inout) :: this
1087 integer :: n_points
1088 !!Perhaps this should be kind dp
1089 !!this is to get around that x,y,z often is 4 dimensional...
1090 !!Should maybe add interface for 1d aswell
1091 real(kind=rp) :: x(n_points,1,1,1)
1092 real(kind=rp) :: y(n_points,1,1,1)
1093 real(kind=rp) :: z(n_points,1,1,1)
1094 integer :: i
1095
1096 call this%free_points()
1097 call this%free_points_local()
1098
1099 this%n_points = n_points
1100
1102
1103 !Deepcopy of coordinates
1104 do i = 1, n_points
1105 this%xyz(1, i) = x(i,1,1,1)
1106 this%xyz(2, i) = y(i,1,1,1)
1107 this%xyz(3, i) = z(i,1,1,1)
1108 end do
1109
1111
1120 subroutine global_interpolation_find_coords1d(this, x, y, z, n_points)
1121 class(global_interpolation_t), intent(inout) :: this
1122 integer :: n_points
1123 real(kind=rp) :: x(n_points)
1124 real(kind=rp) :: y(n_points)
1125 real(kind=rp) :: z(n_points)
1126 integer :: i
1127
1128 call this%free_points()
1129 call this%free_points_local()
1130
1131 this%n_points = n_points
1132
1134
1135 !Deepcopy of coordinates
1136 do i = 1, n_points
1137 this%xyz(1, i) = x(i)
1138 this%xyz(2, i) = y(i)
1139 this%xyz(3, i) = z(i)
1140 end do
1141
1143
1145
1146
1148 class(global_interpolation_t) :: this
1149
1150 if (allocated(this%xyz)) deallocate(this%xyz)
1151 if (allocated(this%rst)) deallocate(this%rst)
1152 if (allocated(this%pe_owner)) deallocate(this%pe_owner)
1153 if (allocated(this%el_owner0)) deallocate(this%el_owner0)
1154
1155 allocate(this%pe_owner(this%n_points))
1156 allocate(this%el_owner0(this%n_points))
1157 allocate(this%xyz(3, this%n_points))
1158 allocate(this%rst(3, this%n_points))
1159 if (neko_bcknd_device .eq. 1) then
1160 call device_map(this%el_owner0, this%el_owner0_d, this%n_points)
1161 end if
1162
1164
1171 subroutine global_interpolation_find_xyz(this, xyz, n_points)
1172 class(global_interpolation_t), intent(inout) :: this
1173 integer, intent(in) :: n_points
1174 !!Perhaps this should be kind dp
1175 real(kind=rp), intent(inout) :: xyz(3, n_points)
1176
1177
1178 call this%free_points()
1179 call this%free_points_local()
1180
1181 this%n_points = n_points
1182
1184
1186 call copy(this%xyz, xyz, 3 * n_points)
1187
1189
1190 end subroutine global_interpolation_find_xyz
1191
1200 subroutine global_interpolation_find_and_redist(this, xyz, n_points)
1201 class(global_interpolation_t), intent(inout) :: this
1202 integer, intent(inout) :: n_points
1203 !!Perhaps this should be kind dp
1204 real(kind=rp), allocatable, intent(inout) :: xyz(:,:)
1205
1206
1207 call this%free_points()
1208 call this%free_points_local()
1209
1210
1211 this%n_points = n_points
1213
1215 call copy(this%xyz, xyz, 3 * n_points)
1216
1218 call this%free_points()
1219 this%n_points = this%n_points_local
1220 n_points = this%n_points_local
1222 if (allocated(xyz)) then
1223 deallocate(xyz)
1224 end if
1225 allocate(xyz(3, n_points))
1226
1227 call copy(xyz, this%xyz_local, 3*n_points)
1228 call copy(this%rst, this%rst_local, 3*n_points)
1229 call copy(this%xyz, this%xyz_local, 3*n_points)
1230 this%pe_owner = this%pe_rank
1231 this%el_owner0 = this%el_owner0_local
1232 if (neko_bcknd_device .eq. 1) then
1233 call device_memcpy(this%el_owner0, this%el_owner0_d, &
1234 this%n_points, host_to_device, sync = .true.)
1235 end if
1236 this%all_points_local = .true.
1237
1239
1246 subroutine global_interpolation_evaluate_masked(this, interp_values, &
1247 field, mask, on_host)
1248 class(global_interpolation_t), target, intent(inout) :: this
1249 real(kind=rp), intent(inout), target :: interp_values(this%n_points)
1250 real(kind=rp), intent(inout), target :: field(this%n_dof)
1251 type(mask_t), intent(in) :: mask
1252 logical, intent(in) :: on_host
1253
1254 call vector_masked_gather_copy(this%masked_field, field, mask, this%n_dof)
1255 call this%evaluate(interp_values, this%masked_field%x, on_host)
1256
1258
1263 subroutine global_interpolation_evaluate(this, interp_values, field, on_host)
1264 class(global_interpolation_t), target, intent(inout) :: this
1265 real(kind=rp), intent(inout), target :: interp_values(this%n_points)
1266 real(kind=rp), intent(inout), target :: field(this%nelv*this%Xh%lxyz)
1267 logical, intent(in) :: on_host
1268 type(c_ptr) :: interp_d
1269
1270 if (.not. this%all_points_local) then
1271 call this%local_interp%evaluate(this%temp_local%x, &
1272 this%el_owner0_local, field, this%nelv, on_host)
1273 if (neko_bcknd_device .eq. 1 .and. .not. on_host) then
1274 call device_memcpy(this%temp_local%x, this%temp_local%x_d, &
1275 this%n_points_local, device_to_host, .true.)
1276 end if
1277 interp_values = 0.0_rp
1278 call this%glb_intrp_comm%sendrecv(this%temp_local%x, interp_values, &
1279 this%n_points_local, this%n_points)
1280 if (neko_bcknd_device .eq. 1 .and. .not. on_host .and. &
1281 this%n_points .gt. 0) then
1282 interp_d = device_get_ptr(interp_values)
1283 call device_memcpy(interp_values, interp_d, &
1284 this%n_points, host_to_device, .false.)
1285 end if
1286 else
1287 call this%local_interp%evaluate(interp_values, this%el_owner0_local, &
1288 field, this%nelv, on_host)
1289 end if
1290
1291 end subroutine global_interpolation_evaluate
1292
1293
1304 function rst_cmp(rst1, rst2, res1, res2, tol) result(rst2_better)
1305 real(kind=rp) :: rst1(3), res1(3)
1306 real(kind=rp) :: rst2(3), res2(3)
1307 real(kind=dp) :: tol
1308 logical :: rst2_better
1309 !If rst1 is invalid and rst2 is valid, take rst2
1310 ! If both invalidl, take smallest residual
1311 if (abs(rst1(1)) .gt. 1.0_xp+tol .or. &
1312 abs(rst1(2)) .gt. 1.0_xp+tol .or. &
1313 abs(rst1(3)) .gt. 1.0_xp+tol) then
1314 if (abs(rst2(1)) .le. 1.0_xp+tol .and. &
1315 abs(rst2(2)) .le. 1.0_xp+tol .and. &
1316 abs(rst2(3)) .le. 1.0_xp+tol) then
1317 rst2_better = .true.
1318 else
1319 rst2_better = (norm2(real(res2, xp)) .lt. norm2(real(res1, xp)))
1320 end if
1321 else
1323 rst2_better = (norm2(real(res2, xp)) .lt. norm2(real(res1, xp)) .and. &
1324 abs(rst2(1)) .le. 1.0_xp+tol .and. &
1325 abs(rst2(2)) .le. 1.0_xp+tol .and. &
1326 abs(rst2(3)) .le. 1.0_xp+tol)
1327 end if
1328 end function rst_cmp
1329
1330end module global_interpolation
double real
Return the device pointer for an associated Fortran array.
Definition device.F90:108
Map a Fortran array to a device (allocate and associate)
Definition device.F90:78
Copy data between host and device (or device and device)
Definition device.F90:72
Unmap a Fortran array from a device (deassociate and free)
Definition device.F90:84
Implements aabb_pe_finder given a dofmap.
integer, parameter, public glob_map_size
Minimum number of total boxes in the aabb tree.
Implements cartesian_pe_finder given a dofmap.
Definition comm.F90:1
integer, public pe_size
MPI size of communicator.
Definition comm.F90:61
integer, public pe_rank
MPI rank.
Definition comm.F90:58
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:45
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:48
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:240
integer, parameter, public device_to_host
Definition device.F90:48
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
Defines global interpolation communication Based on the MPI based gather-scatter kernel.
Implements global_interpolation given a dofmap.
subroutine global_interpolation_free(this)
Destructor.
subroutine global_interpolation_find_coords1d(this, x, y, z, n_points)
Finds the corresponding r,s,t coordinates in the correct global element as well as which process that...
subroutine global_interpolation_init_json_dof(this, dof, params_subdict, comm, mask)
Initialize the global interpolation object on a dofmap.
subroutine global_interpolation_free_points(this)
Destructor for point arrays.
real(kind=dp), parameter, public glob_interp_tol
subroutine global_interpolation_init_point_arrays(this)
subroutine global_interpolation_evaluate(this, interp_values, field, on_host)
Evalute the interpolated value in the points given a field.
subroutine global_interpolation_free_points_local(this)
subroutine global_interpolation_init_dof(this, dof, comm, tol, pad, mask)
Initialize the global interpolation object on a dofmap.
subroutine global_interpolation_find_coords(this, x, y, z, n_points)
Finds the corresponding r,s,t coordinates in the correct global element as well as which process that...
logical function rst_cmp(rst1, rst2, res1, res2, tol)
Compares two sets of rst coordinates and checks whether rst2 is better than rst1 given a tolerance re...
real(kind=dp), parameter, public glob_interp_pad
subroutine global_interpolation_init_xyz(this, x, y, z, gdim, nelv, xh, comm, tol, pad)
Initialize the global interpolation object on a set of coordinates.
subroutine global_interpolation_find_xyz(this, xyz, n_points)
Finds the corresponding r,s,t coordinates in the correct global element as well as which process that...
subroutine global_interpolation_init_json_xyz(this, x, y, z, gdim, nelv, xh, params_subdict, comm)
Initialize the global interpolation object on a set of coordinates, with configuration parameters giv...
subroutine global_interpolation_find_and_redist(this, xyz, n_points)
Finds the corresponding r,s,t coordinates and redistributes the points to the owning rank in the corr...
subroutine global_interpolation_find_common(this)
Common routine for finding the points.
subroutine global_interpolation_check_points(this, x, y, z)
Check the points for validity This is used to check that the points are valid and that the interpolat...
subroutine global_interpolation_evaluate_masked(this, interp_values, field, mask, on_host)
Evaluate the interpolated value in a masked field.
Utilities for retrieving parameters from the case files.
Routines to obtain interpolated values on a set of points with known rst coordinates in elements loca...
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
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:70
Defines a matrix.
Definition matrix.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public xp
Definition num_types.f90:14
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements pe_finder given a dofmap.
Definition pe_finder.f90:35
Defines a function space.
Definition space.f90:34
Implements a dynamic stack ADT.
Definition stack.f90:49
Defines structs that are used... Dont know if we should keep it though.
Definition structs.f90:2
Utilities.
Definition utils.f90:35
subroutine, public vector_masked_gather_copy(a, b, mask, n)
Gather a vector to reduced contigous array .
Defines a vector.
Definition vector.f90:34
Implements global interpolation for arbitrary points in the domain.
Implements global interpolation for arbitrary points in the domain.
Implements global interpolation for arbitrary points in the domain.
Base type for element finder providing element candidates for a given point in the domain.
Definition el_finder.f90:43
Global interpolation communication method.
Implements the settings helper data container for global interpolation.
Implements global interpolation for arbitrary points in the domain.
Type to compute local element (rst) coordinates for a gives set points in physical (xyz) space on a S...
Interpolation on a set of points with known rst coordinates in elements local to this process....
Type for consistently handling masks in Neko. This type encapsulates the mask array and its associate...
Definition mask.f90:51
Implements global interpolation for arbitrary points in the domain.
Definition pe_finder.f90:44
The function space for the SEM solution fields.
Definition space.f90:63
Integer based stack.
Definition stack.f90:77
Pointer to array.
Definition structs.f90:14
#define max(a, b)
Definition tensor.cu:40