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
42 use utils, only: neko_error
52 use el_finder, only: el_finder_t
53 use pe_finder, only: pe_finder_t
54 use comm, only: neko_comm
55 use mpi_f08, only: mpi_sum, mpi_comm, mpi_comm_rank, &
56 mpi_comm_size, mpi_wtime, mpi_allreduce, mpi_in_place, mpi_integer, &
57 mpi_min, mpi_barrier, mpi_reduce_scatter_block, mpi_alltoall, &
58 mpi_isend, mpi_irecv
60 use vector, only: vector_t
62 use matrix, only: matrix_t
63 use math, only: copy, neko_eps
64 use mask, only: mask_t
65 use structs, only : array_ptr_t
66 use, intrinsic :: iso_c_binding, only: c_ptr, c_null_ptr, c_associated
67 implicit none
68 private
69
70 integer, public, parameter :: glob_map_size = 4096
71
73 type, public :: global_interpolation_t
75 type(vector_t) :: x
77 type(vector_t) :: y
79 type(vector_t) :: z
81 integer :: gdim
83 integer :: nelv
85 integer :: glb_nelv
87 type(mpi_comm) :: comm
89 integer :: pe_rank
91 integer :: pe_size
93 type(space_t) :: xh
96 integer :: n_points
98 real(kind=rp), allocatable :: xyz(:,:)
100 real(kind=rp), allocatable :: rst(:,:)
102 integer, allocatable :: pe_owner(:)
104 type(stack_i4_t), allocatable :: points_at_pe(:)
107 integer, allocatable :: el_owner0(:)
108 type(c_ptr) :: el_owner0_d = c_null_ptr
109
111 integer :: n_points_local
112 integer, allocatable :: el_owner0_local(:)
113 type(c_ptr) :: el_owner0_local_d = c_null_ptr
114 real(kind=rp), allocatable :: rst_local(:,:)
115 real(kind=rp), allocatable :: xyz_local(:,:)
117 type(local_interpolator_t) :: local_interp
120 logical :: all_points_local = .false.
122 real(kind=rp) :: tol = neko_eps*1e3_rp
124 real(kind=rp) :: padding = 1e-2_rp
125
129 integer, allocatable :: n_points_pe(:)
130 integer, allocatable :: n_points_offset_pe(:)
133 integer, allocatable :: n_points_pe_local(:)
134 integer, allocatable :: n_points_offset_pe_local(:)
137 class(pe_finder_t), allocatable :: pe_finder
139 class(el_finder_t), allocatable :: el_finder
141 type(legendre_rst_finder_t) :: rst_finder
145 type(vector_t) :: temp_local, temp
146 integer :: n_dof = -1
147 type(vector_t) :: masked_field
148
149 contains
151 procedure, pass(this) :: init_xyz => global_interpolation_init_xyz
153 procedure, pass(this) :: init_dof => global_interpolation_init_dof
155 procedure, pass(this) :: free => global_interpolation_free
157 procedure, pass(this) :: free_points => global_interpolation_free_points
158 procedure, pass(this) :: free_points_local => global_interpolation_free_points_local
159 procedure, pass(this) :: find_points_and_redist => &
164 procedure, pass(this) :: find_points_coords => &
166 procedure, pass(this) :: find_points_coords1d => &
169 procedure, pass(this) :: check_points => &
171 procedure, pass(this) :: find_points_xyz => global_interpolation_find_xyz
172 generic :: find_points => find_points_xyz, find_points_coords, find_points_coords1d
174 procedure, pass(this) :: evaluate => global_interpolation_evaluate
175 procedure, pass(this) :: evaluate_masked => global_interpolation_evaluate_masked
176
178 generic :: init => init_dof, init_xyz
179
181
182contains
183
189 subroutine global_interpolation_init_dof(this, dof, comm, tol, pad, mask)
190 class(global_interpolation_t), target, intent(inout) :: this
191 type(dofmap_t) :: dof
192 type(mpi_comm), optional, intent(in) :: comm
193 real(kind=rp), optional :: tol
194 real(kind=rp), optional :: pad
195 type(mask_t), intent(in), optional :: mask
196 real(kind=rp) :: padding, tolerance
197 integer :: temp_nelv
198
199 if (present(pad)) then
200 padding = pad
201 else
202 padding = 1e-2 ! 1% padding of the bounding boxes
203 end if
204
205 if (present(tol)) then
206 tolerance = tol
207 else
208 tolerance = neko_eps*1e3_xp ! 1% padding of the bounding boxes
209 end if
210
211 ! Store the number of dofs
212 this%n_dof = dof%size()
213 ! NOTE: Passing dof%x(:,1,1,1), etc in init_xyz passes down the entire
214 ! dof%x array and not a slice. It is done this way for
215 ! to get the right dimension (see global_interpolation_init_xyz).
216 if (.not. present(mask)) then
217 call this%init_xyz(dof%x(:,1,1,1), dof%y(:,1,1,1), dof%z(:,1,1,1), &
218 dof%msh%gdim, dof%msh%nelv, dof%Xh, comm,tol = tolerance, pad=padding)
219 else
220
221 ! Initialize a helper field with the size of the mask
222 call this%masked_field%init(mask%size())
223 ! Verify that the mask size is compatible with the dofmap
224 temp_nelv = mask%size() / (dof%Xh%lx*dof%Xh%ly*dof%Xh%lz)
225 if (mod(mask%size(), dof%Xh%lx*dof%Xh%ly*dof%Xh%lz) /= 0) then
226 call neko_error("Mask size must be a multiple of the number of elements in the mesh.")
227 end if
228 ! Initialize with the masked coordinates
229 call this%init_xyz(dof%x(mask%get(),1,1,1), dof%y(mask%get(),1,1,1), dof%z(mask%get(),1,1,1), &
230 dof%msh%gdim, temp_nelv, dof%Xh, comm,tol = tolerance, pad=padding)
231 end if
232
233 end subroutine global_interpolation_init_dof
234
245 subroutine global_interpolation_init_xyz(this, x, y, z, gdim, nelv, Xh, &
246 comm, tol, pad)
247 class(global_interpolation_t), target, intent(inout) :: this
248 real(kind=rp), intent(in) :: x(:)
249 real(kind=rp), intent(in) :: y(:)
250 real(kind=rp), intent(in) :: z(:)
251 integer, intent(in) :: gdim
252 integer, intent(in) :: nelv
253 type(mpi_comm), intent(in), optional :: comm
254 type(space_t), intent(in) :: Xh
255 real(kind=rp), intent(in), optional :: tol
256 real(kind=rp), intent(in), optional :: pad
257 integer :: lx, ly, lz, ierr, i, n
258 character(len=8000) :: log_buf
259 real(kind=dp) :: padding
260 real(kind=rp) :: time1, time_start
261 character(len=255) :: mode_str
262 integer :: boxdim, envvar_len
263
264 call this%free()
265
266 if (present(comm)) then
267 this%comm = comm
268 else
269 this%comm = neko_comm
270 end if
271
272 if (present(pad)) then
273 padding = pad
274 this%padding = pad
275 else
276 padding = 1e-2 ! 1% padding of the bounding boxes
277 this%padding = 1e-2
278 end if
279
280 time_start = mpi_wtime()
281 call mpi_barrier(this%comm)
282
283 call mpi_comm_rank(this%comm, this%pe_rank, ierr)
284 call mpi_comm_size(this%comm, this%pe_size, ierr)
285
286 this%gdim = gdim
287 this%nelv = nelv
288 if (present(tol)) then
289 this%tol = tol
290 else
291 this%tol = neko_eps*1e3_xp
292 end if
293
294 call mpi_allreduce(nelv, this%glb_nelv, 1, mpi_integer, &
295 mpi_sum, this%comm, ierr)
296 lx = xh%lx
297 ly = xh%ly
298 lz = xh%lz
299 n = nelv * lx*ly*lz
300 call this%x%init(n)
301 call this%y%init(n)
302 call this%z%init(n)
303 call copy(this%x%x, x, n)
304 call this%x%copy_from(host_to_device,.false.)
305 call copy(this%y%x, y, n)
306 call this%y%copy_from(host_to_device,.false.)
307 call copy(this%z%x, z, n)
308 call this%z%copy_from(host_to_device,.false.)
309 call this%Xh%init(xh%t, lx, ly, lz)
310
311 ! Initialize n_dof if it was not started with a dof-based constructor
312 if (this%n_dof == -1) then
313 this%n_dof = n
314 end if
315
316 call neko_log%section('Global Interpolation')
317
318 call neko_log%message('Initializing global interpolation')
319 call get_environment_variable("NEKO_GLOBAL_INTERP_EL_FINDER", &
320 mode_str, envvar_len)
321
322 if (envvar_len .gt. 0) then
323 if (mode_str(1:envvar_len) == 'AABB') then
324 allocate(aabb_el_finder_t :: this%el_finder)
325 end if
326 end if
327 call get_environment_variable("NEKO_GLOBAL_INTERP_PE_FINDER", &
328 mode_str, envvar_len)
329
330 if (envvar_len .gt. 0) then
331 if (mode_str(1:envvar_len) == 'AABB') then
332 allocate(aabb_pe_finder_t :: this%pe_finder)
333 end if
334 end if
335
336 if (.not. allocated(this%el_finder)) then
337 allocate(cartesian_el_finder_t :: this%el_finder)
338 end if
339 if (.not. allocated(this%pe_finder)) then
340 allocate(cartesian_pe_finder_t :: this%pe_finder)
341 end if
342 select type(el_find => this%el_finder)
343 type is (aabb_el_finder_t)
344 call neko_log%message('Using AABB element finder')
345 call el_find%init(x, y, z, nelv, xh, padding)
346 type is (cartesian_el_finder_t)
347 call neko_log%message('Using Cartesian element finder')
348 boxdim = max(lx*int(real(nelv,xp)**(1.0_xp/3.0_xp)),2)
349 boxdim = min(boxdim, 300)
350 call el_find%init(x, y, z, nelv, xh, boxdim, padding)
351 class default
352 call neko_error('Unknown element finder type')
353 end select
354
355 select type(pe_find => this%pe_finder)
356 type is (aabb_pe_finder_t)
357 call neko_log%message('Using AABB PE finder')
358 call pe_find%init(this%x%x, this%y%x, this%z%x, &
359 nelv, xh, this%comm, padding)
360 type is (cartesian_pe_finder_t)
361 call neko_log%message('Using Cartesian PE finder')
362 boxdim = lx*int(real(this%glb_nelv,xp)**(1.0_xp/3.0_xp))
363 boxdim = max(boxdim,32)
364 boxdim = min(boxdim, &
365 int(8.0_xp*(30000.0_xp*this%pe_size)**(1.0_xp/3.0_xp)))
366 call pe_find%init(this%x%x, this%y%x, this%z%x, &
367 nelv, xh, this%comm, boxdim, padding)
368 class default
369 call neko_error('Unknown PE finder type')
370 end select
371
372 call this%rst_finder%init(this%x%x, this%y%x, this%z%x, nelv, xh, this%tol)
373 if (allocated(this%n_points_pe)) deallocate(this%n_points_pe)
374 if (allocated(this%n_points_pe_local)) deallocate(this%n_points_pe_local)
375 if (allocated(this%n_points_offset_pe_local)) &
376 deallocate(this%n_points_offset_pe_local)
377 if (allocated(this%n_points_offset_pe)) deallocate(this%n_points_offset_pe)
378 allocate(this%n_points_pe(0:(this%pe_size-1)))
379 allocate(this%n_points_offset_pe(0:(this%pe_size-1)))
380 allocate(this%n_points_pe_local(0:(this%pe_size-1)))
381 allocate(this%n_points_offset_pe_local(0:(this%pe_size-1)))
382 allocate(this%points_at_pe(0:(this%pe_size-1)))
383 do i = 0, this%pe_size-1
384 call this%points_at_pe(i)%init()
385 end do
386 call mpi_barrier(this%comm)
387 time1 = mpi_wtime()
388 write(log_buf, '(A,E15.7)') &
389 'Global interpolation initialized (s):', time1-time_start
390 call neko_log%message(log_buf)
391 call neko_log%end_section()
392 end subroutine global_interpolation_init_xyz
393
394
397 class(global_interpolation_t), target, intent(inout) :: this
398 integer :: i
399
400 call this%x%free()
401 call this%y%free()
402 call this%z%free()
403 call this%Xh%free()
404
405 this%nelv = 0
406 this%gdim = 0
407
408 call this%free_points()
409 call this%free_points_local()
410 call this%local_interp%free()
411 if (allocated(this%el_finder)) then
412 call this%el_finder%free()
413 deallocate(this%el_finder)
414 end if
415 if (allocated(this%pe_finder)) then
416 call this%pe_finder%free()
417 deallocate(this%pe_finder)
418 end if
419 call this%rst_finder%free()
420
421 call this%temp_local%free()
422 call this%temp%free()
423 if (allocated(this%points_at_pe)) then
424 do i = 0, this%pe_size-1
425 call this%points_at_pe(i)%free()
426 end do
427 deallocate(this%points_at_pe)
428 end if
429 if (allocated(this%n_points_pe)) deallocate(this%n_points_pe)
430 if (allocated(this%n_points_pe_local)) deallocate(this%n_points_pe_local)
431 if (allocated(this%n_points_offset_pe_local)) &
432 deallocate(this%n_points_offset_pe_local)
433 if (allocated(this%n_points_offset_pe)) deallocate(this%n_points_offset_pe)
434
435
436
437 end subroutine global_interpolation_free
438
441 class(global_interpolation_t), target, intent(inout) :: this
442
443 this%n_points = 0
444 this%all_points_local = .false.
445
446 if (allocated(this%xyz)) deallocate(this%xyz)
447 if (allocated(this%rst)) deallocate(this%rst)
448 if (allocated(this%pe_owner)) deallocate(this%pe_owner)
449 if (allocated(this%el_owner0)) deallocate(this%el_owner0)
450
451 if (c_associated(this%el_owner0_d)) then
452 call device_free(this%el_owner0_d)
453 end if
454
455 call this%glb_intrp_comm%free()
456
457
459
460
462 class(global_interpolation_t), target, intent(inout) :: this
463
464 this%n_points_local = 0
465 this%all_points_local = .false.
466
467 if (allocated(this%xyz_local)) deallocate(this%xyz_local)
468 if (allocated(this%rst_local)) deallocate(this%rst_local)
469 if (allocated(this%el_owner0_local)) deallocate(this%el_owner0_local)
470
471 if (c_associated(this%el_owner0_local_d)) then
472 call device_free(this%el_owner0_local_d)
473 end if
474
476
477
480 class(global_interpolation_t), target, intent(inout) :: this
481 character(len=8000) :: log_buf
482 type(vector_t) :: x_t
483 type(vector_t) :: y_t
484 type(vector_t) :: z_t
485 type(matrix_t) :: rst_local_cand
486 type(vector_t) :: resx
487 type(vector_t) :: resy
488 type(vector_t) :: resz
489 type(c_ptr) :: el_cands_d
490 type(matrix_t) :: res
491 integer :: i, j, stupid_intent
492 type(stack_i4_t), target :: all_el_candidates
493 integer, allocatable :: n_el_cands(:)
494 integer, contiguous, pointer :: el_cands(:), point_ids(:), send_recv(:)
495 real(kind=rp), allocatable :: res_results(:,:)
496 real(kind=rp), allocatable :: rst_results(:,:)
497 integer, allocatable :: el_owner_results(:)
498 integer :: ierr, ii, n_point_cand, n_glb_point_cand, point_id, rank
499 real(kind=rp) :: time1, time2, time_start
500 !Temp stuff for glb_intrp_comm
501 type(stack_i4_t) :: send_pe, recv_pe
502 type(glb_intrp_comm_t) :: glb_intrp_find, glb_intrp_find_back
503 type(stack_i4_t) :: send_pe_find, recv_pe_find
504
505 el_cands_d = c_null_ptr
506 call neko_log%section('Global Interpolation')
507 call glb_intrp_find%init_dofs(this%pe_size)
508 call send_pe_find%init()
509 call recv_pe_find%init()
510 call mpi_barrier(this%comm)
511 time_start = mpi_wtime()
512 write(log_buf,'(A)') 'Global interpolation, finding points'
513 call neko_log%message(log_buf)
514 ! Find pe candidates that the points i want may be at
515 ! Add number to n_points_pe_local
516 !Working arrays
517 this%n_points_pe = 0
518 call this%pe_finder%find_batch(this%xyz, this%n_points, &
519 this%points_at_pe, this%n_points_pe)
520 call mpi_barrier(this%comm)
521 time1 = mpi_wtime()
522 write(log_buf, '(A,E15.7)') &
523 'Found PE candidates time since start of findpts (s):', time1-time_start
524 call neko_log%message(log_buf)
525
526 !Send number of points I want to candidates
527 ! n_points_local -> how many points might be at this rank
528 ! n_points_pe_local -> how many points local on this rank that other pes might want
529 this%n_points_pe_local = 0
530 this%n_points_local = 0
531 call mpi_reduce_scatter_block(this%n_points_pe, this%n_points_local, &
532 1, mpi_integer, mpi_sum, this%comm, ierr)
533 call mpi_alltoall(this%n_points_pe, 1, mpi_integer,&
534 this%n_points_pe_local, 1, mpi_integer, this%comm, ierr)
535
536 !Set up offset arrays
537 this%n_points_offset_pe_local(0) = 0
538 this%n_points_offset_pe(0) = 0
539 do i = 1, (this%pe_size - 1)
540 this%n_points_offset_pe_local(i) = this%n_points_pe_local(i-1)&
541 + this%n_points_offset_pe_local(i-1)
542 this%n_points_offset_pe(i) = this%n_points_pe(i-1)&
543 + this%n_points_offset_pe(i-1)
544 end do
545 do i = 0, (this%pe_size-1)
546 if (this%n_points_pe(i) .gt. 0) then
547 call send_pe_find%push(i)
548 point_ids => this%points_at_pe(i)%array()
549 do j = 1, this%n_points_pe(i)
550 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+1)
551 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+2)
552 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+3)
553 end do
554 end if
555 if (this%n_points_pe_local(i) .gt. 0) then
556 call recv_pe_find%push(i)
557 do j = 1, this%n_points_pe_local(i)
558 call glb_intrp_find%recv_dof(i)%push(3*(j+this%n_points_offset_pe_local(i)-1)+1)
559 call glb_intrp_find%recv_dof(i)%push(3*(j+this%n_points_offset_pe_local(i)-1)+2)
560 call glb_intrp_find%recv_dof(i)%push(3*(j+this%n_points_offset_pe_local(i)-1)+3)
561 end do
562 end if
563 end do
564
565
566
567 call glb_intrp_find%init(send_pe_find, recv_pe_find, this%comm)
568
569 call glb_intrp_find_back%init_dofs(this%pe_size)
570 ii = 0
571 do i = 0, (this%pe_size-1)
572 send_recv => glb_intrp_find%recv_dof(i)%array()
573 do j = 1, glb_intrp_find%recv_dof(i)%size()
574 call glb_intrp_find_back%send_dof(i)%push(send_recv(j))
575 end do
576 send_recv => glb_intrp_find%send_dof(i)%array()
577 do j = 1, glb_intrp_find%send_dof(i)%size()
578 ii = ii + 1
579 call glb_intrp_find_back%recv_dof(i)%push(ii)
580 end do
581 end do
582
583 call glb_intrp_find_back%init(recv_pe_find, send_pe_find, this%comm)
584
585
586 if (allocated(this%xyz_local)) then
587 deallocate(this%xyz_local)
588 end if
589 allocate(this%xyz_local(3, this%n_points_local))
590 call glb_intrp_find%sendrecv(this%xyz, this%xyz_local, this%n_points*3, &
591 this%n_points_local*3)
592
593 call mpi_barrier(this%comm)
594 time1 = mpi_wtime()
595 write(log_buf, '(A,E15.7)') &
596 'Sent to points to PE candidates, time since start of find_points (s):', time1-time_start
597 call neko_log%message(log_buf)
598 !Okay, now we need to find the rst...
599 call all_el_candidates%init()
600
601 if (allocated(n_el_cands)) then
602 deallocate(n_el_cands)
603 end if
604
605 allocate(n_el_cands(this%n_points_local))
607 call this%el_finder%find_batch(this%xyz_local, this%n_points_local, &
608 all_el_candidates, n_el_cands)
609
610 n_point_cand = all_el_candidates%size()
611 if (n_point_cand .gt. 1e8) then
612 print *,'Warning, many point candidates on rank', this%pe_rank,'cands:', n_point_cand, &
613 'Consider increasing number of ranks'
614 end if
615 call x_t%init(n_point_cand)
616 call y_t%init(n_point_cand)
617 call z_t%init(n_point_cand)
618 ii = 0
620 do i = 1 , this%n_points_local
621 do j = 1, n_el_cands(i)
622 ii = ii + 1
623 x_t%x(ii) = this%xyz_local(1,i)
624 y_t%x(ii) = this%xyz_local(2,i)
625 z_t%x(ii) = this%xyz_local(3,i)
626 end do
627 end do
628
629 call mpi_barrier(this%comm)
630 time1 = mpi_wtime()
631 write(log_buf, '(A,E15.7)') &
632 'Element candidates found, now time for finding rst,time since start of find_points (s):', time1-time_start
633 call neko_log%message(log_buf)
634 call rst_local_cand%init(3,n_point_cand)
635 call resx%init(n_point_cand)
636 call resy%init(n_point_cand)
637 call resz%init(n_point_cand)
638
639 ! Find rst within all element candidates for target xyz (x_t, y_t, z_t)
640 call mpi_barrier(this%comm)
641 time1 = mpi_wtime()
642 el_cands => all_el_candidates%array()
643 if (neko_bcknd_device .eq. 1) then
644 call x_t%copy_from(host_to_device,.false.)
645
646 call y_t%copy_from(host_to_device,.false.)
647
648 call z_t%copy_from(host_to_device,.false.)
649 call device_map(el_cands, el_cands_d,n_point_cand)
650 call device_memcpy(el_cands, el_cands_d,n_point_cand, &
651 host_to_device, .true.)
652 end if
653
654 call this%rst_finder%find(rst_local_cand, &
655 x_t, y_t, z_t, &
656 el_cands, n_point_cand, &
657 resx, resy, resz)
658 if (neko_bcknd_device .eq. 1) then
659 call rst_local_cand%copy_from(device_to_host,.false.)
660 call resx%copy_from(device_to_host,.false.)
661 call resy%copy_from(device_to_host,.false.)
662 call resz%copy_from(device_to_host,.true.)
663 call device_deassociate(el_cands)
664 call device_free(el_cands_d)
665 end if
666 call mpi_barrier(this%comm)
667
668 time2 = mpi_wtime()
669 write(log_buf, '(A,E15.7)') &
670 'Found rst with Newton iteration, time (s):', time2-time1
671 call neko_log%message(log_buf)
672
673 write(log_buf,'(A,E15.7)') &
674 'Tolerance: ', this%tol
675 call neko_log%message(log_buf)
676 write(log_buf,'(A)') &
677 'Checking validity of points and choosing best candidates.'
678 call neko_log%message(log_buf)
679 call mpi_barrier(this%comm,ierr)
680
681 if (allocated(this%rst_local)) deallocate(this%rst_local)
682 if (allocated(this%el_owner0_local)) deallocate(this%el_owner0_local)
683 allocate(this%rst_local(3,this%n_points_local))
684 allocate(this%el_owner0_local(this%n_points_local))
685 ! Choose the best candidate at this rank
686 ii = 0
687 do i = 1 , this%n_points_local
688 this%xyz_local(1,i) = 10.0
689 this%xyz_local(2,i) = 10.0
690 this%xyz_local(3,i) = 10.0
691 this%rst_local(1,i) = 10.0
692 this%rst_local(2,i) = 10.0
693 this%rst_local(3,i) = 10.0
694 this%el_owner0_local(i) = -1
695 do j = 1, n_el_cands(i)
696 ii = ii + 1
697 if (rst_cmp(this%rst_local(:,i), rst_local_cand%x(:,ii),&
698 this%xyz_local(:,i), (/resx%x(ii),resy%x(ii),resz%x(ii)/), this%padding)) then
699 this%rst_local(1,i) = rst_local_cand%x(1,ii)
700 this%rst_local(2,i) = rst_local_cand%x(2,ii)
701
702 this%rst_local(3,i) = rst_local_cand%x(3,ii)
703 this%xyz_local(1,i) = resx%x(ii)
704 this%xyz_local(2,i) = resy%x(ii)
705 this%xyz_local(3,i) = resz%x(ii)
706 this%el_owner0_local(i) = el_cands(ii)
707 end if
708 ! if (this%pe_rank .eq. 0) print *,i, this%rst_local(:,i), &
709 ! this%xyz_local(:,i), this%el_owner0_local(i)
710 end do
711 end do
712 call res%init(3,this%n_points)
713 n_glb_point_cand = sum(this%n_points_pe)
714 if (allocated(rst_results)) deallocate(rst_results)
715 if (allocated(res_results)) deallocate(res_results)
716 if (allocated(el_owner_results)) deallocate(el_owner_results)
717 allocate(rst_results(3,n_glb_point_cand))
718 allocate(res_results(3,n_glb_point_cand))
719 allocate(el_owner_results(n_glb_point_cand))
720 res = 1e2_rp
721 this%rst = 1e2
722 this%pe_owner = -1
723 this%el_owner0 = -1
724 call glb_intrp_find_back%sendrecv(this%xyz_local, res_results, &
725 this%n_points_local*3, n_glb_point_cand*3)
726 call glb_intrp_find_back%sendrecv(this%rst_local, rst_results, &
727 this%n_points_local*3, n_glb_point_cand*3)
728 do i = 1, size(glb_intrp_find_back%send_pe)
729 rank = glb_intrp_find_back%send_pe(i)
730 call mpi_isend(this%el_owner0_local(this%n_points_offset_pe_local(rank)+1),&
731 this%n_points_pe_local(rank), &
732 mpi_integer, rank, 0, &
733 this%comm, glb_intrp_find_back%send_buf(i)%request, ierr)
734 glb_intrp_find_back%send_buf(i)%flag = .false.
735 end do
736 do i = 1, size(glb_intrp_find_back%recv_pe)
737 rank = glb_intrp_find_back%recv_pe(i)
738 call mpi_irecv(el_owner_results(this%n_points_offset_pe(rank)+1),&
739 this%n_points_pe(rank), &
740 mpi_integer, rank, 0, &
741 this%comm, glb_intrp_find_back%recv_buf(i)%request, ierr)
742 glb_intrp_find_back%recv_buf(i)%flag = .false.
743 end do
744 call glb_intrp_find_back%nbwait_no_op()
745 ii = 0
746 do i = 1, size(glb_intrp_find_back%recv_pe)
747 point_ids => this%points_at_pe(glb_intrp_find_back%recv_pe(i))%array()
748 do j = 1, this%n_points_pe(glb_intrp_find_back%recv_pe(i))
749 point_id = point_ids(j)
750 ii = ii + 1
751 if (rst_cmp(this%rst(:,point_id), rst_results(:,ii), &
752 res%x(:,point_id), res_results(:,ii), this%padding) .or. &
753 this%pe_owner(point_ids(j)) .eq. -1 ) then
754 this%rst(:,point_ids(j)) = rst_results(:,ii)
755 res%x(:,point_ids(j)) = res_results(:,ii)
756 this%pe_owner(point_ids(j)) = glb_intrp_find_back%recv_pe(i)
757 this%el_owner0(point_ids(j)) = el_owner_results(ii)
758 end if
759 ! if (this%pe_rank .eq. 0) print *,point_id, &
760 !this%rst(:,point_ids(j)),res%x(:,point_ids(j)), this%el_owner0(point_ids(j))
761 end do
762 end do
763
764 !OK, now I know the correct rst values of the points I want We now
765 !send the correct rsts to the correct rank (so a point only
766 !belongs to one rank)
767 do i = 0, this%pe_size-1
768 call this%points_at_pe(i)%clear()
769 this%n_points_pe(i) = 0
770 end do
771
772 do i = 1, this%n_points
773 stupid_intent = i
774 if (this%pe_owner(i) .eq. -1 .or. this%el_owner0(i) .eq. -1) then
775 print *, 'No owning rank found for',&
776 ' point ', stupid_intent, ' with coords', this%xyz(:,i), &
777 ' Interpolation will always yield 0.0. Try increase padding.'
778 else
779 call this%points_at_pe(this%pe_owner(i))%push(stupid_intent)
780 this%n_points_pe(this%pe_owner(i)) = this%n_points_pe(this%pe_owner(i)) + 1
781 end if
782 end do
783 call mpi_reduce_scatter_block(this%n_points_pe, this%n_points_local, 1, mpi_integer, &
784 mpi_sum, this%comm, ierr)
785 call mpi_alltoall(this%n_points_pe, 1, mpi_integer,&
786 this%n_points_pe_local, 1, mpi_integer, this%comm, ierr)
787 this%n_points_offset_pe_local(0) = 0
788 this%n_points_offset_pe(0) = 0
789 do i = 1, (this%pe_size - 1)
790 this%n_points_offset_pe_local(i) = this%n_points_pe_local(i-1)&
791 + this%n_points_offset_pe_local(i-1)
792 this%n_points_offset_pe(i) = this%n_points_pe(i-1)&
793 + this%n_points_offset_pe(i-1)
794 end do
795 call send_pe_find%free()
796 call recv_pe_find%free()
797 call glb_intrp_find%free()
798 call send_pe_find%init()
799 call recv_pe_find%init()
800 call glb_intrp_find%init_dofs(this%pe_size)
801 !setup comm to send xyz and rst to chosen ranks
802 do i = 0, (this%pe_size-1)
803 if (this%n_points_pe(i) .gt. 0) then
804 call send_pe_find%push(i)
805 point_ids => this%points_at_pe(i)%array()
806 do j = 1, this%n_points_pe(i)
807 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+1)
808 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+2)
809 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+3)
810 end do
811 end if
812 if (this%n_points_pe_local(i) .gt. 0) then
813 call recv_pe_find%push(i)
814 do j = 1, this%n_points_pe_local(i)
815 call glb_intrp_find%recv_dof(i)%push(3*(j+this%n_points_offset_pe_local(i)-1)+1)
816 call glb_intrp_find%recv_dof(i)%push(3*(j+this%n_points_offset_pe_local(i)-1)+2)
817 call glb_intrp_find%recv_dof(i)%push(3*(j+this%n_points_offset_pe_local(i)-1)+3)
818 end do
819 end if
820 end do
821
822
823 call glb_intrp_find%init(send_pe_find, recv_pe_find, this%comm)
824 call glb_intrp_find%sendrecv(this%xyz, this%xyz_local, this%n_points*3, &
825 this%n_points_local*3)
826 call glb_intrp_find%sendrecv(this%rst, this%rst_local, this%n_points*3, &
827 this%n_points_local*3)
828 ii = 0
829 do i = 1, size(glb_intrp_find%send_pe)
830 rank = glb_intrp_find%send_pe(i)
831 point_ids => this%points_at_pe(rank)%array()
832 do j = 1, this%n_points_pe(rank)
833 ii = ii + 1
834 el_owner_results(ii) = this%el_owner0(point_ids(j))
835 end do
836 call mpi_isend(el_owner_results(this%n_points_offset_pe(rank)+1),&
837 this%n_points_pe(rank), &
838 mpi_integer, rank, 0, &
839 this%comm, glb_intrp_find%send_buf(i)%request, ierr)
840 glb_intrp_find%send_buf(i)%flag = .false.
841 end do
842 do i = 1, size(glb_intrp_find%recv_pe)
843 rank = glb_intrp_find%recv_pe(i)
844 call mpi_irecv(this%el_owner0_local(this%n_points_offset_pe_local(rank)+1), &
845 this%n_points_pe_local(rank), &
846 mpi_integer, rank, 0, &
847 this%comm, glb_intrp_find%recv_buf(i)%request, ierr)
848 glb_intrp_find%recv_buf(i)%flag = .false.
849 end do
850 call glb_intrp_find%nbwait_no_op()
851
852 call glb_intrp_find%free()
853
854 !Set up final way of doing communication
855 call send_pe%init()
856 call recv_pe%init()
857 call this%glb_intrp_comm%init_dofs(this%pe_size)
858 do i = 0, (this%pe_size-1)
859 if (this%n_points_pe(i) .gt. 0) then
860 call recv_pe%push(i)
861 point_ids => this%points_at_pe(i)%array()
862 do j = 1, this%n_points_pe(i)
863 call this%glb_intrp_comm%recv_dof(i)%push(point_ids(j))
864 end do
865 end if
866 if (this%n_points_pe_local(i) .gt. 0) then
867 call send_pe%push(i)
868 do j = 1, this%n_points_pe_local(i)
869 call this%glb_intrp_comm%send_dof(i)%push(j+this%n_points_offset_pe_local(i))
870 end do
871 end if
872 end do
873 call this%glb_intrp_comm%init(send_pe, recv_pe,this%comm)
874
875 !Initialize working arrays for evaluation
876 call this%temp_local%init(this%n_points_local)
877 call this%temp%init(this%n_points)
878
879 !Initialize interpolator for local interpolation
880 call this%local_interp%init(this%Xh, this%rst_local, &
881 this%n_points_local)
882
883
884 if (neko_bcknd_device .eq. 1) then
885 call device_memcpy(this%el_owner0, this%el_owner0_d, &
886 this%n_points, host_to_device, sync = .true.)
887 call device_map(this%el_owner0_local, this%el_owner0_local_d, this%n_points_local)
888 call device_memcpy(this%el_owner0_local, this%el_owner0_local_d, &
889 this%n_points_local, host_to_device, sync = .true.)
890 end if
891
892 call this%check_points(this%x%x,this%y%x, this%z%x)
893
894 !Free stuff
895 call send_pe%free()
896 call recv_pe%free()
897 call glb_intrp_find_back%free()
898 call send_pe_find%free()
899 call recv_pe_find%free()
900 call x_t%free()
901 call y_t%free()
902 call z_t%free()
903 call rst_local_cand%free()
904 call resx%free()
905 call resy%free()
906 call resz%free()
907 call res%free()
908 call all_el_candidates%free()
909
910 if (allocated(n_el_cands)) deallocate(n_el_cands)
911 if (allocated(rst_results)) deallocate(rst_results)
912 if (allocated(res_results)) deallocate(res_results)
913 if (allocated(el_owner_results)) deallocate(el_owner_results)
914 call mpi_barrier(this%comm, ierr)
915 time2 = mpi_wtime()
916 write(log_buf, '(A,E15.7)') 'Global interpolation find points done, time (s):', &
917 time2-time_start
918 call neko_log%message(log_buf)
919 call neko_log%end_section()
920 call neko_log%newline()
922
926 subroutine global_interpolation_check_points(this, x, y, z)
927 class(global_interpolation_t), target, intent(inout) :: this
928 real(kind=rp), intent(inout) :: x(:)
929 real(kind=rp), intent(inout) :: y(:)
930 real(kind=rp), intent(inout) :: z(:)
931 integer :: i, j
932 character(len=8000) :: log_buf
933 real(kind=rp) :: xdiff, ydiff, zdiff
934 logical :: isdiff
935 type(vector_t) :: x_check, y_check, z_check
936
937 call x_check%init(this%n_points)
938 call y_check%init(this%n_points)
939 call z_check%init(this%n_points)
940 call this%evaluate(x_check%x, x, on_host=.true.)
941 call this%evaluate(y_check%x, y, on_host=.true.)
942 call this%evaluate(z_check%x, z, on_host=.true.)
943 write(log_buf,'(A)') 'Checking validity of points.'
944 call neko_log%message(log_buf)
945 j = 0
946 do i = 1 , this%n_points
947 ! Check validity of points
948 isdiff = .false.
949 xdiff = x_check%x(i)-this%xyz(1,i)
950 ydiff = y_check%x(i)-this%xyz(2,i)
951 zdiff = z_check%x(i)-this%xyz(3,i)
952 isdiff = norm2(real((/xdiff,ydiff,zdiff/),xp)) > this%tol
953 if (isdiff) then
954 write(*,*) 'Point ', i,'at rank ', this%pe_rank, 'with coordinates: ', &
955 this%xyz(1, i), this%xyz(2, i), this%xyz(3, i), &
956 'Differ from interpolated coords: ', &
957 x_check%x(i), y_check%x(i), z_check%x(i), &
958 'Actual difference: ', &
959 xdiff, ydiff, zdiff, norm2(real((/xdiff,ydiff,zdiff/),xp)),&
960 'Process, element: ', &
961 this%pe_owner(i), this%el_owner0(i)+1, &
962 'Calculated rst: ', &
963 this%rst(1,i), this%rst(2,i), this%rst(3,i)
964 j = j + 1
965 end if
966 end do
967 call x_check%free()
968 call y_check%free()
969 call z_check%free()
970
972
980 subroutine global_interpolation_find_coords(this, x, y, z, n_points)
981 class(global_interpolation_t), intent(inout) :: this
982 integer :: n_points
983 !!Perhaps this should be kind dp
984 !!this is to get around that x,y,z often is 4 dimensional...
985 !!Should maybe add interface for 1d aswell
986 real(kind=rp) :: x(n_points,1,1,1)
987 real(kind=rp) :: y(n_points,1,1,1)
988 real(kind=rp) :: z(n_points,1,1,1)
989 integer :: i
990
991 call this%free_points()
992 call this%free_points_local()
993
994 this%n_points = n_points
995
997
998 !Deepcopy of coordinates
999 do i = 1, n_points
1000 this%xyz(1, i) = x(i,1,1,1)
1001 this%xyz(2, i) = y(i,1,1,1)
1002 this%xyz(3, i) = z(i,1,1,1)
1003 end do
1004
1006
1015 subroutine global_interpolation_find_coords1d(this, x, y, z, n_points)
1016 class(global_interpolation_t), intent(inout) :: this
1017 integer :: n_points
1018 real(kind=rp) :: x(n_points)
1019 real(kind=rp) :: y(n_points)
1020 real(kind=rp) :: z(n_points)
1021 integer :: i
1022
1023 call this%free_points()
1024 call this%free_points_local()
1025
1026 this%n_points = n_points
1027
1029
1030 !Deepcopy of coordinates
1031 do i = 1, n_points
1032 this%xyz(1, i) = x(i)
1033 this%xyz(2, i) = y(i)
1034 this%xyz(3, i) = z(i)
1035 end do
1036
1038
1040
1041
1043 class(global_interpolation_t) :: this
1044
1045 if (allocated(this%xyz)) deallocate(this%xyz)
1046 if (allocated(this%rst)) deallocate(this%rst)
1047 if (allocated(this%pe_owner)) deallocate(this%pe_owner)
1048 if (allocated(this%el_owner0)) deallocate(this%el_owner0)
1049
1050 allocate(this%pe_owner(this%n_points))
1051 allocate(this%el_owner0(this%n_points))
1052 allocate(this%xyz(3, this%n_points))
1053 allocate(this%rst(3, this%n_points))
1054 if (neko_bcknd_device .eq. 1) then
1055 call device_map(this%el_owner0, this%el_owner0_d, this%n_points)
1056 end if
1057
1059
1066 subroutine global_interpolation_find_xyz(this, xyz, n_points)
1067 class(global_interpolation_t), intent(inout) :: this
1068 integer, intent(in) :: n_points
1069 !!Perhaps this should be kind dp
1070 real(kind=rp), intent(inout) :: xyz(3, n_points)
1071
1072
1073 call this%free_points()
1074 call this%free_points_local()
1075
1076 this%n_points = n_points
1077
1079
1081 call copy(this%xyz, xyz, 3 * n_points)
1082
1084
1085 end subroutine global_interpolation_find_xyz
1086
1095 subroutine global_interpolation_find_and_redist(this, xyz, n_points)
1096 class(global_interpolation_t), intent(inout) :: this
1097 integer, intent(inout) :: n_points
1098 !!Perhaps this should be kind dp
1099 real(kind=rp), allocatable, intent(inout) :: xyz(:,:)
1100
1101
1102 call this%free_points()
1103 call this%free_points_local()
1104
1105
1106 this%n_points = n_points
1108
1110 call copy(this%xyz, xyz, 3 * n_points)
1111
1113 call this%free_points()
1114 this%n_points = this%n_points_local
1115 n_points = this%n_points_local
1117 if (allocated(xyz)) then
1118 deallocate(xyz)
1119 end if
1120 allocate(xyz(3,n_points))
1121
1122 call copy(xyz, this%xyz_local, 3*n_points)
1123 call copy(this%rst, this%rst_local, 3*n_points)
1124 call copy(this%xyz, this%xyz_local, 3*n_points)
1125 this%pe_owner = this%pe_rank
1126 this%el_owner0 = this%el_owner0_local
1127 if (neko_bcknd_device .eq. 1) then
1128 call device_memcpy(this%el_owner0, this%el_owner0_d, &
1129 this%n_points, host_to_device, sync = .true.)
1130 end if
1131 this%all_points_local = .true.
1132
1134
1141 subroutine global_interpolation_evaluate_masked(this, interp_values, &
1142 field, mask, on_host)
1143 class(global_interpolation_t), target, intent(inout) :: this
1144 real(kind=rp), intent(inout), target :: interp_values(this%n_points)
1145 real(kind=rp), intent(inout), target :: field(this%n_dof)
1146 type(mask_t), intent(in) :: mask
1147 logical, intent(in) :: on_host
1148
1149 call vector_masked_gather_copy(this%masked_field, field, mask, this%n_dof)
1150 call this%evaluate(interp_values, this%masked_field%x, on_host)
1151
1153
1158 subroutine global_interpolation_evaluate(this, interp_values, field, on_host)
1159 class(global_interpolation_t), target, intent(inout) :: this
1160 real(kind=rp), intent(inout), target :: interp_values(this%n_points)
1161 real(kind=rp), intent(inout), target :: field(this%nelv*this%Xh%lxyz)
1162 logical, intent(in) :: on_host
1163 type(c_ptr) :: interp_d
1164
1165 if (.not. this%all_points_local) then
1166 call this%local_interp%evaluate(this%temp_local%x, &
1167 this%el_owner0_local, field, this%nelv, on_host)
1168 if (neko_bcknd_device .eq. 1 .and. .not. on_host) then
1169 call device_memcpy(this%temp_local%x, this%temp_local%x_d, &
1170 this%n_points_local, device_to_host, .true.)
1171 end if
1172 interp_values = 0.0_rp
1173 call this%glb_intrp_comm%sendrecv(this%temp_local%x, interp_values, &
1174 this%n_points_local, this%n_points)
1175 if (neko_bcknd_device .eq. 1 .and. .not. on_host) then
1176 interp_d = device_get_ptr(interp_values)
1177 call device_memcpy(interp_values, interp_d, &
1178 this%n_points, host_to_device, .false.)
1179 end if
1180 else
1181 call this%local_interp%evaluate(interp_values, this%el_owner0_local, &
1182 field, this%nelv, on_host)
1183 end if
1184
1185 end subroutine global_interpolation_evaluate
1186
1187
1198 function rst_cmp(rst1, rst2,res1, res2, tol) result(rst2_better)
1199 real(kind=rp) :: rst1(3), res1(3)
1200 real(kind=rp) :: rst2(3), res2(3)
1201 real(kind=rp) :: tol
1202 logical :: rst2_better
1203 !If rst1 is invalid and rst2 is valid, take rst2
1204 ! If both invalidl, take smallest residual
1205 if (abs(rst1(1)) .gt. 1.0_xp+tol .or. &
1206 abs(rst1(2)) .gt. 1.0_xp+tol .or. &
1207 abs(rst1(3)) .gt. 1.0_xp+tol) then
1208 if (abs(rst2(1)) .le. 1.0_xp+tol .and. &
1209 abs(rst2(2)) .le. 1.0_xp+tol .and. &
1210 abs(rst2(3)) .le. 1.0_xp+tol) then
1211 rst2_better = .true.
1212 else
1213 rst2_better = (norm2(real(res2,xp)) .lt. norm2(real(res1,xp)))
1214 end if
1215 else
1217 rst2_better = (norm2(real(res2,xp)) .lt. norm2(real(res1,xp)) .and.&
1218 abs(rst2(1)) .le. 1.0_xp+tol .and. &
1219 abs(rst2(2)) .le. 1.0_xp+tol .and. &
1220 abs(rst2(3)) .le. 1.0_xp+tol)
1221 end if
1222 end function rst_cmp
1223
1224end module global_interpolation
double real
Deassociate a Fortran array from a device pointer.
Definition device.F90:95
Return the device pointer for an associated Fortran array.
Definition device.F90:101
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:219
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_free_points(this)
Destructor for point arrays.
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...
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_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.
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:249
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:69
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:35
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 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:63
Pointer to array.
Definition structs.f90:14
#define max(a, b)
Definition tensor.cu:40