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)) deallocate(this%el_owner0)
537
538 if (c_associated(this%el_owner0_d)) then
539 call device_free(this%el_owner0_d)
540 end if
541
542 call this%glb_intrp_comm%free()
543
544
546
547
549 class(global_interpolation_t), target, intent(inout) :: this
550
551 this%n_points_local = 0
552 this%all_points_local = .false.
553
554 if (allocated(this%xyz_local)) deallocate(this%xyz_local)
555 if (allocated(this%rst_local)) deallocate(this%rst_local)
556 if (allocated(this%el_owner0_local)) deallocate(this%el_owner0_local)
557
558 if (c_associated(this%el_owner0_local_d)) then
559 call device_free(this%el_owner0_local_d)
560 end if
561
563
564
567 class(global_interpolation_t), target, intent(inout) :: this
568 character(len=8000) :: log_buf
569 type(vector_t) :: x_t
570 type(vector_t) :: y_t
571 type(vector_t) :: z_t
572 type(matrix_t) :: rst_local_cand
573 type(vector_t) :: resx
574 type(vector_t) :: resy
575 type(vector_t) :: resz
576 type(c_ptr) :: el_cands_d
577 type(matrix_t) :: res
578 integer :: i, j, stupid_intent
579 type(stack_i4_t), target :: all_el_candidates
580 integer, allocatable :: n_el_cands(:)
581 integer, contiguous, pointer :: el_cands(:), point_ids(:), send_recv(:)
582 real(kind=rp), allocatable :: res_results(:,:)
583 real(kind=rp), allocatable :: rst_results(:,:)
584 integer, allocatable :: el_owner_results(:)
585 integer :: ierr, ii, n_point_cand, n_glb_point_cand, point_id, rank
586 real(kind=rp) :: time1, time2, time_start
587 !Temp stuff for glb_intrp_comm
588 type(stack_i4_t) :: send_pe, recv_pe
589 type(glb_intrp_comm_t) :: glb_intrp_find, glb_intrp_find_back
590 type(stack_i4_t) :: send_pe_find, recv_pe_find
591
592 el_cands_d = c_null_ptr
593 call neko_log%section('Global Interpolation')
594 call glb_intrp_find%init_dofs(this%pe_size)
595 call send_pe_find%init()
596 call recv_pe_find%init()
597 call mpi_barrier(this%comm)
598 time_start = mpi_wtime()
599 write(log_buf, '(A)') 'Global interpolation, finding points'
600 call neko_log%message(log_buf)
601 ! Find pe candidates that the points i want may be at
602 ! Add number to n_points_pe_local
603 !Working arrays
604 this%n_points_pe = 0
605 call this%pe_finder%find_batch(this%xyz, this%n_points, &
606 this%points_at_pe, this%n_points_pe)
607 call mpi_barrier(this%comm)
608 time1 = mpi_wtime()
609 write(log_buf, '(A,E15.7)') &
610 'Found PE candidates time since start of findpts (s):', &
611 time1 - time_start
612 call neko_log%message(log_buf)
613
614 !Send number of points I want to candidates
615 ! n_points_local -> how many points might be at this rank
616 ! n_points_pe_local -> how many points local on this rank that other pes
617 ! might want
618 this%n_points_pe_local = 0
619 this%n_points_local = 0
620 call mpi_reduce_scatter_block(this%n_points_pe, this%n_points_local, &
621 1, mpi_integer, mpi_sum, this%comm, ierr)
622 call mpi_alltoall(this%n_points_pe, 1, mpi_integer,&
623 this%n_points_pe_local, 1, mpi_integer, this%comm, ierr)
624
625 !Set up offset arrays
626 this%n_points_offset_pe_local(0) = 0
627 this%n_points_offset_pe(0) = 0
628 do i = 1, (this%pe_size - 1)
629 this%n_points_offset_pe_local(i) = this%n_points_pe_local(i-1)&
630 + this%n_points_offset_pe_local(i-1)
631 this%n_points_offset_pe(i) = this%n_points_pe(i-1)&
632 + this%n_points_offset_pe(i-1)
633 end do
634 do i = 0, (this%pe_size-1)
635 if (this%n_points_pe(i) .gt. 0) then
636 call send_pe_find%push(i)
637 point_ids => this%points_at_pe(i)%array()
638 do j = 1, this%n_points_pe(i)
639 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+1)
640 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+2)
641 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+3)
642 end do
643 end if
644 if (this%n_points_pe_local(i) .gt. 0) then
645 call recv_pe_find%push(i)
646 do j = 1, this%n_points_pe_local(i)
647 call glb_intrp_find%recv_dof(i)%push(3*(j + &
648 this%n_points_offset_pe_local(i) - 1) + 1)
649 call glb_intrp_find%recv_dof(i)%push(3*(j + &
650 this%n_points_offset_pe_local(i) - 1) + 2)
651 call glb_intrp_find%recv_dof(i)%push(3*(j + &
652 this%n_points_offset_pe_local(i) - 1) + 3)
653 end do
654 end if
655 end do
656
657
658
659 call glb_intrp_find%init(send_pe_find, recv_pe_find, this%comm)
660
661 call glb_intrp_find_back%init_dofs(this%pe_size)
662 ii = 0
663 do i = 0, (this%pe_size-1)
664 send_recv => glb_intrp_find%recv_dof(i)%array()
665 do j = 1, glb_intrp_find%recv_dof(i)%size()
666 call glb_intrp_find_back%send_dof(i)%push(send_recv(j))
667 end do
668 send_recv => glb_intrp_find%send_dof(i)%array()
669 do j = 1, glb_intrp_find%send_dof(i)%size()
670 ii = ii + 1
671 call glb_intrp_find_back%recv_dof(i)%push(ii)
672 end do
673 end do
674
675 call glb_intrp_find_back%init(recv_pe_find, send_pe_find, this%comm)
676
677
678 if (allocated(this%xyz_local)) then
679 deallocate(this%xyz_local)
680 end if
681 allocate(this%xyz_local(3, this%n_points_local))
682 call glb_intrp_find%sendrecv(this%xyz, this%xyz_local, this%n_points*3, &
683 this%n_points_local*3)
684
685 call mpi_barrier(this%comm)
686 time1 = mpi_wtime()
687 write(log_buf, '(A,E15.7)') &
688 'Sent to points to PE candidates, time since start of ' &
689 // 'find_points (s):', time1 - time_start
690 call neko_log%message(log_buf)
691
692 !Okay, now we need to find the rst...
693 call all_el_candidates%init()
694
695 if (allocated(n_el_cands)) then
696 deallocate(n_el_cands)
697 end if
698
699 allocate(n_el_cands(this%n_points_local))
701 call this%el_finder%find_batch(this%xyz_local, this%n_points_local, &
702 all_el_candidates, n_el_cands)
703
704 n_point_cand = all_el_candidates%size()
705 if (n_point_cand .gt. 1e8) then
706 print *,'Warning, many point candidates on rank', this%pe_rank, &
707 'cands:', n_point_cand, &
708 'Consider increasing number of ranks'
709 end if
710 call x_t%init(n_point_cand)
711 call y_t%init(n_point_cand)
712 call z_t%init(n_point_cand)
713 ii = 0
715 do i = 1 , this%n_points_local
716 do j = 1, n_el_cands(i)
717 ii = ii + 1
718 x_t%x(ii) = this%xyz_local(1,i)
719 y_t%x(ii) = this%xyz_local(2,i)
720 z_t%x(ii) = this%xyz_local(3,i)
721 end do
722 end do
723
724 call mpi_barrier(this%comm)
725 time1 = mpi_wtime()
726 write(log_buf, '(A,E15.7)') &
727 'Element candidates found, now time for finding rst,time ' // &
728 'since start of find_points (s):', time1 - time_start
729 call neko_log%message(log_buf)
730 call rst_local_cand%init(3,n_point_cand)
731 call resx%init(n_point_cand)
732 call resy%init(n_point_cand)
733 call resz%init(n_point_cand)
734
735 ! Find rst within all element candidates for target xyz (x_t, y_t, z_t)
736 call mpi_barrier(this%comm)
737 time1 = mpi_wtime()
738 el_cands => all_el_candidates%array()
739 if (neko_bcknd_device .eq. 1) then
740 call x_t%copy_from(host_to_device,.false.)
741
742 call y_t%copy_from(host_to_device,.false.)
743
744 call z_t%copy_from(host_to_device,.false.)
745 call device_map(el_cands, el_cands_d,n_point_cand)
746 call device_memcpy(el_cands, el_cands_d,n_point_cand, &
747 host_to_device, .true.)
748 end if
749
750 call this%rst_finder%find(rst_local_cand, &
751 x_t, y_t, z_t, &
752 el_cands, n_point_cand, &
753 resx, resy, resz)
754 if (neko_bcknd_device .eq. 1) then
755 call rst_local_cand%copy_from(device_to_host,.false.)
756 call resx%copy_from(device_to_host,.false.)
757 call resy%copy_from(device_to_host,.false.)
758 call resz%copy_from(device_to_host,.true.)
759 call device_deassociate(el_cands)
760 call device_free(el_cands_d)
761 end if
762 call mpi_barrier(this%comm)
763
764 time2 = mpi_wtime()
765 write(log_buf, '(A,E15.7)') &
766 'Found rst with Newton iteration, time (s):', time2-time1
767 call neko_log%message(log_buf)
768
769 write(log_buf, '(A)') &
770 'Checking validity of points and choosing best candidates.'
771 call neko_log%message(log_buf)
772 call mpi_barrier(this%comm,ierr)
773
774 if (allocated(this%rst_local)) deallocate(this%rst_local)
775 if (allocated(this%el_owner0_local)) deallocate(this%el_owner0_local)
776 allocate(this%rst_local(3,this%n_points_local))
777 allocate(this%el_owner0_local(this%n_points_local))
778 ! Choose the best candidate at this rank
779 ii = 0
780 do i = 1 , this%n_points_local
781 this%xyz_local(1,i) = 10.0
782 this%xyz_local(2,i) = 10.0
783 this%xyz_local(3,i) = 10.0
784 this%rst_local(1,i) = 10.0
785 this%rst_local(2,i) = 10.0
786 this%rst_local(3,i) = 10.0
787 this%el_owner0_local(i) = -1
788 do j = 1, n_el_cands(i)
789 ii = ii + 1
790 if (rst_cmp(this%rst_local(:,i), rst_local_cand%x(:,ii),&
791 this%xyz_local(:,i), [resx%x(ii),resy%x(ii),resz%x(ii)], &
792 this%padding)) then
793 this%rst_local(1,i) = rst_local_cand%x(1,ii)
794 this%rst_local(2,i) = rst_local_cand%x(2,ii)
795
796 this%rst_local(3,i) = rst_local_cand%x(3,ii)
797 this%xyz_local(1,i) = resx%x(ii)
798 this%xyz_local(2,i) = resy%x(ii)
799 this%xyz_local(3,i) = resz%x(ii)
800 this%el_owner0_local(i) = el_cands(ii)
801 end if
802 ! if (this%pe_rank .eq. 0) print *,i, this%rst_local(:,i), &
803 ! this%xyz_local(:,i), this%el_owner0_local(i)
804 end do
805 end do
806 call res%init(3,this%n_points)
807 n_glb_point_cand = sum(this%n_points_pe)
808 if (allocated(rst_results)) deallocate(rst_results)
809 if (allocated(res_results)) deallocate(res_results)
810 if (allocated(el_owner_results)) deallocate(el_owner_results)
811 allocate(rst_results(3,n_glb_point_cand))
812 allocate(res_results(3,n_glb_point_cand))
813 allocate(el_owner_results(n_glb_point_cand))
814 res = 1e2_rp
815 this%rst = 1e2
816 this%pe_owner = -1
817 this%el_owner0 = -1
818 call glb_intrp_find_back%sendrecv(this%xyz_local, res_results, &
819 this%n_points_local*3, n_glb_point_cand*3)
820 call glb_intrp_find_back%sendrecv(this%rst_local, rst_results, &
821 this%n_points_local*3, n_glb_point_cand*3)
822 do i = 1, size(glb_intrp_find_back%send_pe)
823 rank = glb_intrp_find_back%send_pe(i)
824 call mpi_isend(this%el_owner0_local( &
825 this%n_points_offset_pe_local(rank) + 1), &
826 this%n_points_pe_local(rank), &
827 mpi_integer, rank, 0, &
828 this%comm, glb_intrp_find_back%send_buf(i)%request, ierr)
829 glb_intrp_find_back%send_buf(i)%flag = .false.
830 end do
831 do i = 1, size(glb_intrp_find_back%recv_pe)
832 rank = glb_intrp_find_back%recv_pe(i)
833 call mpi_irecv(el_owner_results(this%n_points_offset_pe(rank)+1),&
834 this%n_points_pe(rank), &
835 mpi_integer, rank, 0, &
836 this%comm, glb_intrp_find_back%recv_buf(i)%request, ierr)
837 glb_intrp_find_back%recv_buf(i)%flag = .false.
838 end do
839 call glb_intrp_find_back%nbwait_no_op()
840 ii = 0
841 do i = 1, size(glb_intrp_find_back%recv_pe)
842 point_ids => this%points_at_pe(glb_intrp_find_back%recv_pe(i))%array()
843 do j = 1, this%n_points_pe(glb_intrp_find_back%recv_pe(i))
844 point_id = point_ids(j)
845 ii = ii + 1
846 if (rst_cmp(this%rst(:,point_id), rst_results(:,ii), &
847 res%x(:,point_id), res_results(:,ii), this%padding) .or. &
848 this%pe_owner(point_ids(j)) .eq. -1 ) then
849 this%rst(:,point_ids(j)) = rst_results(:,ii)
850 res%x(:,point_ids(j)) = res_results(:,ii)
851 this%pe_owner(point_ids(j)) = glb_intrp_find_back%recv_pe(i)
852 this%el_owner0(point_ids(j)) = el_owner_results(ii)
853 end if
854 ! if (this%pe_rank .eq. 0) print *,point_id, &
855 !this%rst(:,point_ids(j)),res%x(:,point_ids(j)),
856 ! this%el_owner0(point_ids(j))
857 end do
858 end do
859
860 !OK, now I know the correct rst values of the points I want We now
861 !send the correct rsts to the correct rank (so a point only
862 !belongs to one rank)
863 do i = 0, this%pe_size-1
864 call this%points_at_pe(i)%clear()
865 this%n_points_pe(i) = 0
866 end do
867
868 do i = 1, this%n_points
869 stupid_intent = i
870 if (this%pe_owner(i) .eq. -1 .or. this%el_owner0(i) .eq. -1) then
871 print *, 'No owning rank found for',&
872 ' point ', stupid_intent, ' with coords', this%xyz(:,i), &
873 ' Interpolation will always yield 0.0. Try increase padding.'
874 else
875 call this%points_at_pe(this%pe_owner(i))%push(stupid_intent)
876 this%n_points_pe(this%pe_owner(i)) = &
877 this%n_points_pe(this%pe_owner(i)) + 1
878 end if
879 end do
880 call mpi_reduce_scatter_block(this%n_points_pe, this%n_points_local, 1, &
881 mpi_integer, mpi_sum, this%comm, ierr)
882 call mpi_alltoall(this%n_points_pe, 1, mpi_integer, &
883 this%n_points_pe_local, 1, mpi_integer, this%comm, ierr)
884 this%n_points_offset_pe_local(0) = 0
885 this%n_points_offset_pe(0) = 0
886 do i = 1, (this%pe_size - 1)
887 this%n_points_offset_pe_local(i) = this%n_points_pe_local(i-1)&
888 + this%n_points_offset_pe_local(i-1)
889 this%n_points_offset_pe(i) = this%n_points_pe(i-1)&
890 + this%n_points_offset_pe(i-1)
891 end do
892 call send_pe_find%free()
893 call recv_pe_find%free()
894 call glb_intrp_find%free()
895 call send_pe_find%init()
896 call recv_pe_find%init()
897 call glb_intrp_find%init_dofs(this%pe_size)
898 !setup comm to send xyz and rst to chosen ranks
899 do i = 0, (this%pe_size-1)
900 if (this%n_points_pe(i) .gt. 0) then
901 call send_pe_find%push(i)
902 point_ids => this%points_at_pe(i)%array()
903 do j = 1, this%n_points_pe(i)
904 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 1)
905 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 2)
906 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 3)
907 end do
908 end if
909 if (this%n_points_pe_local(i) .gt. 0) then
910 call recv_pe_find%push(i)
911 do j = 1, this%n_points_pe_local(i)
912 call glb_intrp_find%recv_dof(i)%push(3*(j + &
913 this%n_points_offset_pe_local(i) - 1) + 1)
914 call glb_intrp_find%recv_dof(i)%push(3*(j + &
915 this%n_points_offset_pe_local(i) - 1) + 2)
916 call glb_intrp_find%recv_dof(i)%push(3*(j + &
917 this%n_points_offset_pe_local(i) - 1) + 3)
918 end do
919 end if
920 end do
921
922
923 call glb_intrp_find%init(send_pe_find, recv_pe_find, this%comm)
924 call glb_intrp_find%sendrecv(this%xyz, this%xyz_local, this%n_points*3, &
925 this%n_points_local*3)
926 call glb_intrp_find%sendrecv(this%rst, this%rst_local, this%n_points*3, &
927 this%n_points_local*3)
928 ii = 0
929 do i = 1, size(glb_intrp_find%send_pe)
930 rank = glb_intrp_find%send_pe(i)
931 point_ids => this%points_at_pe(rank)%array()
932 do j = 1, this%n_points_pe(rank)
933 ii = ii + 1
934 el_owner_results(ii) = this%el_owner0(point_ids(j))
935 end do
936 call mpi_isend(el_owner_results(this%n_points_offset_pe(rank) + 1),&
937 this%n_points_pe(rank), &
938 mpi_integer, rank, 0, &
939 this%comm, glb_intrp_find%send_buf(i)%request, ierr)
940 glb_intrp_find%send_buf(i)%flag = .false.
941 end do
942 do i = 1, size(glb_intrp_find%recv_pe)
943 rank = glb_intrp_find%recv_pe(i)
944 call mpi_irecv(this%el_owner0_local( &
945 this%n_points_offset_pe_local(rank) + 1), &
946 this%n_points_pe_local(rank), &
947 mpi_integer, rank, 0, &
948 this%comm, glb_intrp_find%recv_buf(i)%request, ierr)
949 glb_intrp_find%recv_buf(i)%flag = .false.
950 end do
951 call glb_intrp_find%nbwait_no_op()
952
953 call glb_intrp_find%free()
954
955 !Set up final way of doing communication
956 call send_pe%init()
957 call recv_pe%init()
958 call this%glb_intrp_comm%init_dofs(this%pe_size)
959 do i = 0, (this%pe_size-1)
960 if (this%n_points_pe(i) .gt. 0) then
961 call recv_pe%push(i)
962 point_ids => this%points_at_pe(i)%array()
963 do j = 1, this%n_points_pe(i)
964 call this%glb_intrp_comm%recv_dof(i)%push(point_ids(j))
965 end do
966 end if
967 if (this%n_points_pe_local(i) .gt. 0) then
968 call send_pe%push(i)
969 do j = 1, this%n_points_pe_local(i)
970 call this%glb_intrp_comm%send_dof(i)%push(j + &
971 this%n_points_offset_pe_local(i))
972 end do
973 end if
974 end do
975 call this%glb_intrp_comm%init(send_pe, recv_pe,this%comm)
976
977 !Initialize working arrays for evaluation
978 call this%temp_local%init(this%n_points_local)
979 call this%temp%init(this%n_points)
980
981 !Initialize interpolator for local interpolation
982 call this%local_interp%init(this%Xh, this%rst_local, &
983 this%n_points_local)
984
985
986 if (neko_bcknd_device .eq. 1) then
987 call device_memcpy(this%el_owner0, this%el_owner0_d, &
988 this%n_points, host_to_device, sync = .true.)
989 call device_map(this%el_owner0_local, this%el_owner0_local_d, &
990 this%n_points_local)
991 call device_memcpy(this%el_owner0_local, this%el_owner0_local_d, &
992 this%n_points_local, host_to_device, sync = .true.)
993 end if
994
995 call this%check_points(this%x%x,this%y%x, this%z%x)
996
997 !Free stuff
998 call send_pe%free()
999 call recv_pe%free()
1000 call glb_intrp_find_back%free()
1001 call send_pe_find%free()
1002 call recv_pe_find%free()
1003 call x_t%free()
1004 call y_t%free()
1005 call z_t%free()
1006 call rst_local_cand%free()
1007 call resx%free()
1008 call resy%free()
1009 call resz%free()
1010 call res%free()
1011 call all_el_candidates%free()
1012
1013 if (allocated(n_el_cands)) deallocate(n_el_cands)
1014 if (allocated(rst_results)) deallocate(rst_results)
1015 if (allocated(res_results)) deallocate(res_results)
1016 if (allocated(el_owner_results)) deallocate(el_owner_results)
1017 call mpi_barrier(this%comm, ierr)
1018 time2 = mpi_wtime()
1019 write(log_buf, '(A,E15.7)') 'Global interpolation find points ' // &
1020 'done, time (s):', time2-time_start
1021 call neko_log%message(log_buf)
1022 call neko_log%end_section()
1023 call neko_log%newline()
1025
1029 subroutine global_interpolation_check_points(this, x, y, z)
1030 class(global_interpolation_t), target, intent(inout) :: this
1031 real(kind=rp), intent(inout) :: x(:)
1032 real(kind=rp), intent(inout) :: y(:)
1033 real(kind=rp), intent(inout) :: z(:)
1034 integer :: i, j
1035 character(len=8000) :: log_buf
1036 real(kind=rp) :: xdiff, ydiff, zdiff
1037 logical :: isdiff
1038 type(vector_t) :: x_check, y_check, z_check
1039
1040 call x_check%init(this%n_points)
1041 call y_check%init(this%n_points)
1042 call z_check%init(this%n_points)
1043 call this%evaluate(x_check%x, x, on_host = .true.)
1044 call this%evaluate(y_check%x, y, on_host = .true.)
1045 call this%evaluate(z_check%x, z, on_host = .true.)
1046 write(log_buf, '(A)') 'Checking validity of points.'
1047 call neko_log%message(log_buf)
1048 j = 0
1049 do i = 1 , this%n_points
1050 ! Check validity of points
1051 isdiff = .false.
1052 xdiff = x_check%x(i)-this%xyz(1,i)
1053 ydiff = y_check%x(i)-this%xyz(2,i)
1054 zdiff = z_check%x(i)-this%xyz(3,i)
1055 isdiff = norm2(real([xdiff,ydiff,zdiff],xp)) > this%tolerance
1056 if (isdiff) then
1057 write(*,*) 'Point ', i,'at rank ', this%pe_rank, &
1058 'with coordinates: ', &
1059 this%xyz(1, i), this%xyz(2, i), this%xyz(3, i), &
1060 'Differ from interpolated coords: ', &
1061 x_check%x(i), y_check%x(i), z_check%x(i), &
1062 'Actual difference: ', &
1063 xdiff, ydiff, zdiff, norm2(real([xdiff,ydiff,zdiff],xp)),&
1064 'Process, element: ', &
1065 this%pe_owner(i), this%el_owner0(i)+1, &
1066 'Calculated rst: ', &
1067 this%rst(1,i), this%rst(2,i), this%rst(3,i)
1068 j = j + 1
1069 end if
1070 end do
1071 call x_check%free()
1072 call y_check%free()
1073 call z_check%free()
1074
1076
1084 subroutine global_interpolation_find_coords(this, x, y, z, n_points)
1085 class(global_interpolation_t), intent(inout) :: this
1086 integer :: n_points
1087 !!Perhaps this should be kind dp
1088 !!this is to get around that x,y,z often is 4 dimensional...
1089 !!Should maybe add interface for 1d aswell
1090 real(kind=rp) :: x(n_points,1,1,1)
1091 real(kind=rp) :: y(n_points,1,1,1)
1092 real(kind=rp) :: z(n_points,1,1,1)
1093 integer :: i
1094
1095 call this%free_points()
1096 call this%free_points_local()
1097
1098 this%n_points = n_points
1099
1101
1102 !Deepcopy of coordinates
1103 do i = 1, n_points
1104 this%xyz(1, i) = x(i,1,1,1)
1105 this%xyz(2, i) = y(i,1,1,1)
1106 this%xyz(3, i) = z(i,1,1,1)
1107 end do
1108
1110
1119 subroutine global_interpolation_find_coords1d(this, x, y, z, n_points)
1120 class(global_interpolation_t), intent(inout) :: this
1121 integer :: n_points
1122 real(kind=rp) :: x(n_points)
1123 real(kind=rp) :: y(n_points)
1124 real(kind=rp) :: z(n_points)
1125 integer :: i
1126
1127 call this%free_points()
1128 call this%free_points_local()
1129
1130 this%n_points = n_points
1131
1133
1134 !Deepcopy of coordinates
1135 do i = 1, n_points
1136 this%xyz(1, i) = x(i)
1137 this%xyz(2, i) = y(i)
1138 this%xyz(3, i) = z(i)
1139 end do
1140
1142
1144
1145
1147 class(global_interpolation_t) :: this
1148
1149 if (allocated(this%xyz)) deallocate(this%xyz)
1150 if (allocated(this%rst)) deallocate(this%rst)
1151 if (allocated(this%pe_owner)) deallocate(this%pe_owner)
1152 if (allocated(this%el_owner0)) deallocate(this%el_owner0)
1153
1154 allocate(this%pe_owner(this%n_points))
1155 allocate(this%el_owner0(this%n_points))
1156 allocate(this%xyz(3, this%n_points))
1157 allocate(this%rst(3, this%n_points))
1158 if (neko_bcknd_device .eq. 1) then
1159 call device_map(this%el_owner0, this%el_owner0_d, this%n_points)
1160 end if
1161
1163
1170 subroutine global_interpolation_find_xyz(this, xyz, n_points)
1171 class(global_interpolation_t), intent(inout) :: this
1172 integer, intent(in) :: n_points
1173 !!Perhaps this should be kind dp
1174 real(kind=rp), intent(inout) :: xyz(3, n_points)
1175
1176
1177 call this%free_points()
1178 call this%free_points_local()
1179
1180 this%n_points = n_points
1181
1183
1185 call copy(this%xyz, xyz, 3 * n_points)
1186
1188
1189 end subroutine global_interpolation_find_xyz
1190
1199 subroutine global_interpolation_find_and_redist(this, xyz, n_points)
1200 class(global_interpolation_t), intent(inout) :: this
1201 integer, intent(inout) :: n_points
1202 !!Perhaps this should be kind dp
1203 real(kind=rp), allocatable, intent(inout) :: xyz(:,:)
1204
1205
1206 call this%free_points()
1207 call this%free_points_local()
1208
1209
1210 this%n_points = n_points
1212
1214 call copy(this%xyz, xyz, 3 * n_points)
1215
1217 call this%free_points()
1218 this%n_points = this%n_points_local
1219 n_points = this%n_points_local
1221 if (allocated(xyz)) then
1222 deallocate(xyz)
1223 end if
1224 allocate(xyz(3,n_points))
1225
1226 call copy(xyz, this%xyz_local, 3*n_points)
1227 call copy(this%rst, this%rst_local, 3*n_points)
1228 call copy(this%xyz, this%xyz_local, 3*n_points)
1229 this%pe_owner = this%pe_rank
1230 this%el_owner0 = this%el_owner0_local
1231 if (neko_bcknd_device .eq. 1) then
1232 call device_memcpy(this%el_owner0, this%el_owner0_d, &
1233 this%n_points, host_to_device, sync = .true.)
1234 end if
1235 this%all_points_local = .true.
1236
1238
1245 subroutine global_interpolation_evaluate_masked(this, interp_values, &
1246 field, mask, on_host)
1247 class(global_interpolation_t), target, intent(inout) :: this
1248 real(kind=rp), intent(inout), target :: interp_values(this%n_points)
1249 real(kind=rp), intent(inout), target :: field(this%n_dof)
1250 type(mask_t), intent(in) :: mask
1251 logical, intent(in) :: on_host
1252
1253 call vector_masked_gather_copy(this%masked_field, field, mask, this%n_dof)
1254 call this%evaluate(interp_values, this%masked_field%x, on_host)
1255
1257
1262 subroutine global_interpolation_evaluate(this, interp_values, field, on_host)
1263 class(global_interpolation_t), target, intent(inout) :: this
1264 real(kind=rp), intent(inout), target :: interp_values(this%n_points)
1265 real(kind=rp), intent(inout), target :: field(this%nelv*this%Xh%lxyz)
1266 logical, intent(in) :: on_host
1267 type(c_ptr) :: interp_d
1268
1269 if (.not. this%all_points_local) then
1270 call this%local_interp%evaluate(this%temp_local%x, &
1271 this%el_owner0_local, field, this%nelv, on_host)
1272 if (neko_bcknd_device .eq. 1 .and. .not. on_host) then
1273 call device_memcpy(this%temp_local%x, this%temp_local%x_d, &
1274 this%n_points_local, device_to_host, .true.)
1275 end if
1276 interp_values = 0.0_rp
1277 call this%glb_intrp_comm%sendrecv(this%temp_local%x, interp_values, &
1278 this%n_points_local, this%n_points)
1279 if (neko_bcknd_device .eq. 1 .and. .not. on_host .and. this%n_points .gt. 0) then
1280 interp_d = device_get_ptr(interp_values)
1281 call device_memcpy(interp_values, interp_d, &
1282 this%n_points, host_to_device, .false.)
1283 end if
1284 else
1285 call this%local_interp%evaluate(interp_values, this%el_owner0_local, &
1286 field, this%nelv, on_host)
1287 end if
1288
1289 end subroutine global_interpolation_evaluate
1290
1291
1302 function rst_cmp(rst1, rst2,res1, res2, tol) result(rst2_better)
1303 real(kind=rp) :: rst1(3), res1(3)
1304 real(kind=rp) :: rst2(3), res2(3)
1305 real(kind=dp) :: tol
1306 logical :: rst2_better
1307 !If rst1 is invalid and rst2 is valid, take rst2
1308 ! If both invalidl, take smallest residual
1309 if (abs(rst1(1)) .gt. 1.0_xp+tol .or. &
1310 abs(rst1(2)) .gt. 1.0_xp+tol .or. &
1311 abs(rst1(3)) .gt. 1.0_xp+tol) then
1312 if (abs(rst2(1)) .le. 1.0_xp+tol .and. &
1313 abs(rst2(2)) .le. 1.0_xp+tol .and. &
1314 abs(rst2(3)) .le. 1.0_xp+tol) then
1315 rst2_better = .true.
1316 else
1317 rst2_better = (norm2(real(res2,xp)) .lt. norm2(real(res1,xp)))
1318 end if
1319 else
1321 rst2_better = (norm2(real(res2,xp)) .lt. norm2(real(res1,xp)) .and.&
1322 abs(rst2(1)) .le. 1.0_xp+tol .and. &
1323 abs(rst2(2)) .le. 1.0_xp+tol .and. &
1324 abs(rst2(3)) .le. 1.0_xp+tol)
1325 end if
1326 end function rst_cmp
1327
1328end module global_interpolation
double real
Deassociate a Fortran array from a device pointer.
Definition device.F90:101
Return the device pointer for an associated Fortran array.
Definition device.F90:107
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
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:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:225
integer, parameter, public device_to_host
Definition device.F90:47
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:79
integer, parameter, public log_size
Definition log.f90:45
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:252
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