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
77 type, public :: global_interpolation_t
79 type(vector_t) :: x
81 type(vector_t) :: y
83 type(vector_t) :: z
85 integer :: gdim
87 integer :: nelv
89 integer :: glb_nelv
91 type(mpi_comm) :: comm
93 integer :: pe_rank
95 integer :: pe_size
97 type(space_t) :: xh
100 integer :: n_points
102 real(kind=rp), allocatable :: xyz(:,:)
104 real(kind=rp), allocatable :: rst(:,:)
106 integer, allocatable :: pe_owner(:)
108 type(stack_i4_t), allocatable :: points_at_pe(:)
111 integer, allocatable :: el_owner0(:)
112 type(c_ptr) :: el_owner0_d = c_null_ptr
113
115 integer :: n_points_local
116 integer, allocatable :: el_owner0_local(:)
117 type(c_ptr) :: el_owner0_local_d = c_null_ptr
118 real(kind=rp), allocatable :: rst_local(:,:)
119 real(kind=rp), allocatable :: xyz_local(:,:)
121 type(local_interpolator_t) :: local_interp
124 logical :: all_points_local = .false.
126 real(kind=dp) :: tolerance = glob_interp_tol
128 real(kind=dp) :: padding = glob_interp_pad
129
133 integer, allocatable :: n_points_pe(:)
134 integer, allocatable :: n_points_offset_pe(:)
137 integer, allocatable :: n_points_pe_local(:)
138 integer, allocatable :: n_points_offset_pe_local(:)
141 class(pe_finder_t), allocatable :: pe_finder
143 class(el_finder_t), allocatable :: el_finder
145 type(legendre_rst_finder_t) :: rst_finder
150 type(vector_t) :: temp_local, temp
151 integer :: n_dof = -1
152 type(vector_t) :: masked_field
153
154 contains
157 procedure, pass(this) :: init_json_xyz => &
161 procedure, pass(this) :: init_json_dof => &
165 procedure, pass(this) :: init_xyz => global_interpolation_init_xyz
167 procedure, pass(this) :: init_dof => global_interpolation_init_dof
169 procedure, pass(this) :: free => global_interpolation_free
171 procedure, pass(this) :: free_points => global_interpolation_free_points
172 procedure, pass(this) :: free_points_local => &
174 procedure, pass(this) :: find_points_and_redist => &
179 procedure, pass(this) :: find_points_coords => &
181 procedure, pass(this) :: find_points_coords1d => &
184 procedure, pass(this) :: check_points => &
186 procedure, pass(this) :: find_points_xyz => global_interpolation_find_xyz
187 generic :: find_points => find_points_xyz, find_points_coords, &
188 find_points_coords1d
190 procedure, pass(this) :: evaluate => global_interpolation_evaluate
191 procedure, pass(this) :: evaluate_masked => &
193
195 generic :: init => init_dof, init_xyz, init_json_xyz, init_json_dof
196
198
199contains
200
214 subroutine global_interpolation_init_json_xyz(this, x, y, z, gdim, nelv, Xh, &
215 params_subdict, comm)
216 class(global_interpolation_t), target, intent(inout) :: this
217 real(kind=rp), intent(in) :: x(:)
218 real(kind=rp), intent(in) :: y(:)
219 real(kind=rp), intent(in) :: z(:)
220 integer, intent(in) :: gdim
221 integer, intent(in) :: nelv
222 type(space_t), intent(in) :: Xh
223 type(json_file), intent(inout) :: params_subdict
224 type(mpi_comm), intent(in), optional :: comm
225
226 real(kind=dp) :: tol, pad
227
228 call json_get_or_lookup_or_default(params_subdict, 'tolerance', &
229 tol, glob_interp_tol)
230 call json_get_or_lookup_or_default(params_subdict, 'padding', &
231 pad, glob_interp_pad)
232
233 call this%init_xyz(x, y, z, gdim, nelv, xh, comm = comm, tol = tol, &
234 pad = pad)
235
237
245 subroutine global_interpolation_init_json_dof(this, dof, params_subdict, &
246 comm, mask)
247 class(global_interpolation_t), target, intent(inout) :: this
248 type(dofmap_t) :: dof
249 type(json_file), intent(inout) :: params_subdict
250 type(mpi_comm), optional, intent(in) :: comm
251 type(mask_t), intent(in), optional :: mask
252
253 real(kind=dp) :: tol, pad
254
255 call json_get_or_lookup_or_default(params_subdict, 'tolerance', &
256 tol, glob_interp_tol)
257 call json_get_or_lookup_or_default(params_subdict, 'padding', &
258 pad, glob_interp_pad)
259
260 call this%init_dof(dof, comm = comm, tol = tol, pad = pad, mask = mask)
261
263
270 subroutine global_interpolation_init_dof(this, dof, comm, tol, pad, mask)
271 class(global_interpolation_t), target, intent(inout) :: this
272 type(dofmap_t) :: dof
273 type(mpi_comm), optional, intent(in) :: comm
274 real(kind=dp), optional :: tol
275 real(kind=dp), optional :: pad
276 type(mask_t), intent(in), optional :: mask
277
278 integer :: temp_nelv
279
280 ! Store the number of dofs
281 this%n_dof = dof%size()
282 ! NOTE: Passing dof%x(:,1,1,1), etc in init_xyz passes down the entire
283 ! dof%x array and not a slice. It is done this way for
284 ! to get the right dimension (see global_interpolation_init_xyz).
285 if (.not. present(mask)) then
286 call this%init_xyz(dof%x(:,1,1,1), dof%y(:,1,1,1), dof%z(:,1,1,1), &
287 dof%msh%gdim, dof%msh%nelv, dof%Xh, comm = comm, &
288 tol = tol, pad = pad)
289 else
290
291 ! Initialize a helper field with the size of the mask
292 call this%masked_field%init(mask%size())
293 ! Verify that the mask size is compatible with the dofmap
294 temp_nelv = mask%size() / (dof%Xh%lx*dof%Xh%ly*dof%Xh%lz)
295 if (mod(mask%size(), dof%Xh%lx*dof%Xh%ly*dof%Xh%lz) /= 0) then
296 call neko_error("Mask size must be a multiple of the number of" // &
297 " elements in the mesh.")
298 end if
299 ! Initialize with the masked coordinates
300 call this%init_xyz(dof%x(mask%get(),1,1,1), dof%y(mask%get(),1,1,1), &
301 dof%z(mask%get(),1,1,1), dof%msh%gdim, temp_nelv, dof%Xh, &
302 comm = comm, tol = tol, pad = pad)
303 end if
304
305 end subroutine global_interpolation_init_dof
306
317 subroutine global_interpolation_init_xyz(this, x, y, z, gdim, nelv, Xh, &
318 comm, tol, pad)
319 class(global_interpolation_t), target, intent(inout) :: this
320 real(kind=rp), intent(in) :: x(:)
321 real(kind=rp), intent(in) :: y(:)
322 real(kind=rp), intent(in) :: z(:)
323 integer, intent(in) :: gdim
324 integer, intent(in) :: nelv
325 type(space_t), intent(in) :: Xh
326 type(mpi_comm), intent(in), optional :: comm
327 real(kind=dp), intent(in), optional :: tol
328 real(kind=dp), intent(in), optional :: pad
329
330 integer :: lx, ly, lz, ierr, i, n
331 character(len=8000) :: log_buf
332 !real(kind=dp) :: padding ! <-- note that this padding is in dp
333 real(kind=rp) :: time1, time_start
334 character(len=255) :: mode_str
335 integer :: boxdim, envvar_len
336
337 call neko_log%section('Global Interpolation')
338 call neko_log%message('Initializing global interpolation')
339
340 call this%free()
341
342 ! Set communicator
343 if (present(comm)) then
344 this%comm = comm
345 else
346 this%comm = neko_comm
347 end if
348
349 ! Set point search parameters
350 this%padding = glob_interp_pad
351 if (present(pad)) this%padding = pad
352
353 this%tolerance = glob_interp_tol
354 if (present(tol)) this%tolerance = tol
355
356 write(log_buf, '(A,E15.7)') &
357 'Tolerance: ', this%tolerance
358 call neko_log%message(log_buf)
359 write(log_buf, '(A,E15.7)') &
360 'Padding : ', this%padding
361 call neko_log%message(log_buf)
362
363 time_start = mpi_wtime()
364 call mpi_barrier(this%comm)
365
366 call mpi_comm_rank(this%comm, this%pe_rank, ierr)
367 call mpi_comm_size(this%comm, this%pe_size, ierr)
368
369 this%gdim = gdim
370 this%nelv = nelv
371
372 call mpi_allreduce(nelv, this%glb_nelv, 1, mpi_integer, &
373 mpi_sum, this%comm, ierr)
374 lx = xh%lx
375 ly = xh%ly
376 lz = xh%lz
377 n = nelv * lx*ly*lz
378 call this%x%init(n)
379 call this%y%init(n)
380 call this%z%init(n)
381 call copy(this%x%x, x, n)
382 call this%x%copy_from(host_to_device,.false.)
383 call copy(this%y%x, y, n)
384 call this%y%copy_from(host_to_device,.false.)
385 call copy(this%z%x, z, n)
386 call this%z%copy_from(host_to_device,.false.)
387 call this%Xh%init(xh%t, lx, ly, lz)
388
389 ! Initialize n_dof if it was not started with a dof-based constructor
390 if (this%n_dof == -1) then
391 this%n_dof = n
392 end if
393
394
395 call get_environment_variable("NEKO_GLOBAL_INTERP_EL_FINDER", &
396 mode_str, envvar_len)
397
398 if (envvar_len .gt. 0) then
399 if (mode_str(1:envvar_len) == 'AABB') then
400 allocate(aabb_el_finder_t :: this%el_finder)
401 end if
402 end if
403 call get_environment_variable("NEKO_GLOBAL_INTERP_PE_FINDER", &
404 mode_str, envvar_len)
405
406 if (envvar_len .gt. 0) then
407 if (mode_str(1:envvar_len) == 'AABB') then
408 allocate(aabb_pe_finder_t :: this%pe_finder)
409 end if
410 end if
411
412 if (.not. allocated(this%el_finder)) then
413 allocate(cartesian_el_finder_t :: this%el_finder)
414 end if
415 if (.not. allocated(this%pe_finder)) then
416 allocate(cartesian_pe_finder_t :: this%pe_finder)
417 end if
418 select type (el_find => this%el_finder)
419 type is (aabb_el_finder_t)
420 call neko_log%message('Using AABB element finder')
421 call el_find%init(x, y, z, nelv, xh, this%padding)
422 type is (cartesian_el_finder_t)
423 call neko_log%message('Using Cartesian element finder')
424 boxdim = max(lx*int(real(nelv,xp)**(1.0_xp/3.0_xp)),2)
425 boxdim = min(boxdim, 300)
426 call el_find%init(x, y, z, nelv, xh, boxdim, this%padding)
427 class default
428 call neko_error('Unknown element finder type')
429 end select
430
431 select type (pe_find => this%pe_finder)
432 type is (aabb_pe_finder_t)
433 call neko_log%message('Using AABB PE finder')
434 call pe_find%init(this%x%x, this%y%x, this%z%x, &
435 nelv, xh, this%comm, this%padding)
436 type is (cartesian_pe_finder_t)
437 call neko_log%message('Using Cartesian PE finder')
438 boxdim = lx*int(real(this%glb_nelv,xp)**(1.0_xp/3.0_xp))
439 boxdim = max(boxdim,32)
440 boxdim = min(boxdim, &
441 int(8.0_xp*(30000.0_xp*this%pe_size)**(1.0_xp/3.0_xp)))
442 call pe_find%init(this%x%x, this%y%x, this%z%x, &
443 nelv, xh, this%comm, boxdim, this%padding)
444 class default
445 call neko_error('Unknown PE finder type')
446 end select
447
448 call this%rst_finder%init(this%x%x, this%y%x, this%z%x, nelv, xh, &
449 this%tolerance)
450 if (allocated(this%n_points_pe)) deallocate(this%n_points_pe)
451 if (allocated(this%n_points_pe_local)) deallocate(this%n_points_pe_local)
452 if (allocated(this%n_points_offset_pe_local)) &
453 deallocate(this%n_points_offset_pe_local)
454 if (allocated(this%n_points_offset_pe)) deallocate(this%n_points_offset_pe)
455 allocate(this%n_points_pe(0:(this%pe_size-1)))
456 allocate(this%n_points_offset_pe(0:(this%pe_size-1)))
457 allocate(this%n_points_pe_local(0:(this%pe_size-1)))
458 allocate(this%n_points_offset_pe_local(0:(this%pe_size-1)))
459 allocate(this%points_at_pe(0:(this%pe_size-1)))
460 do i = 0, this%pe_size-1
461 call this%points_at_pe(i)%init()
462 end do
463 call mpi_barrier(this%comm)
464 time1 = mpi_wtime()
465 write(log_buf, '(A,E15.7)') &
466 'Global interpolation initialized (s):', time1-time_start
467 call neko_log%message(log_buf)
468 call neko_log%end_section()
469 end subroutine global_interpolation_init_xyz
470
471
474 class(global_interpolation_t), target, intent(inout) :: this
475 integer :: i
476
477 call this%x%free()
478 call this%y%free()
479 call this%z%free()
480 call this%Xh%free()
481
482 this%nelv = 0
483 this%gdim = 0
484
485 call this%free_points()
486 call this%free_points_local()
487 call this%local_interp%free()
488 if (allocated(this%el_finder)) then
489 call this%el_finder%free()
490 deallocate(this%el_finder)
491 end if
492 if (allocated(this%pe_finder)) then
493 call this%pe_finder%free()
494 deallocate(this%pe_finder)
495 end if
496 call this%rst_finder%free()
497
498 call this%temp_local%free()
499 call this%temp%free()
500 if (allocated(this%points_at_pe)) then
501 do i = 0, this%pe_size-1
502 call this%points_at_pe(i)%free()
503 end do
504 deallocate(this%points_at_pe)
505 end if
506 if (allocated(this%n_points_pe)) deallocate(this%n_points_pe)
507 if (allocated(this%n_points_pe_local)) deallocate(this%n_points_pe_local)
508 if (allocated(this%n_points_offset_pe_local)) &
509 deallocate(this%n_points_offset_pe_local)
510 if (allocated(this%n_points_offset_pe)) deallocate(this%n_points_offset_pe)
511
512
513
514 end subroutine global_interpolation_free
515
518 class(global_interpolation_t), target, intent(inout) :: this
519
520 this%n_points = 0
521 this%all_points_local = .false.
522
523 if (allocated(this%xyz)) deallocate(this%xyz)
524 if (allocated(this%rst)) deallocate(this%rst)
525 if (allocated(this%pe_owner)) deallocate(this%pe_owner)
526 if (allocated(this%el_owner0)) deallocate(this%el_owner0)
527
528 if (c_associated(this%el_owner0_d)) then
529 call device_free(this%el_owner0_d)
530 end if
531
532 call this%glb_intrp_comm%free()
533
534
536
537
539 class(global_interpolation_t), target, intent(inout) :: this
540
541 this%n_points_local = 0
542 this%all_points_local = .false.
543
544 if (allocated(this%xyz_local)) deallocate(this%xyz_local)
545 if (allocated(this%rst_local)) deallocate(this%rst_local)
546 if (allocated(this%el_owner0_local)) deallocate(this%el_owner0_local)
547
548 if (c_associated(this%el_owner0_local_d)) then
549 call device_free(this%el_owner0_local_d)
550 end if
551
553
554
557 class(global_interpolation_t), target, intent(inout) :: this
558 character(len=8000) :: log_buf
559 type(vector_t) :: x_t
560 type(vector_t) :: y_t
561 type(vector_t) :: z_t
562 type(matrix_t) :: rst_local_cand
563 type(vector_t) :: resx
564 type(vector_t) :: resy
565 type(vector_t) :: resz
566 type(c_ptr) :: el_cands_d
567 type(matrix_t) :: res
568 integer :: i, j, stupid_intent
569 type(stack_i4_t), target :: all_el_candidates
570 integer, allocatable :: n_el_cands(:)
571 integer, contiguous, pointer :: el_cands(:), point_ids(:), send_recv(:)
572 real(kind=rp), allocatable :: res_results(:,:)
573 real(kind=rp), allocatable :: rst_results(:,:)
574 integer, allocatable :: el_owner_results(:)
575 integer :: ierr, ii, n_point_cand, n_glb_point_cand, point_id, rank
576 real(kind=rp) :: time1, time2, time_start
577 !Temp stuff for glb_intrp_comm
578 type(stack_i4_t) :: send_pe, recv_pe
579 type(glb_intrp_comm_t) :: glb_intrp_find, glb_intrp_find_back
580 type(stack_i4_t) :: send_pe_find, recv_pe_find
581
582 el_cands_d = c_null_ptr
583 call neko_log%section('Global Interpolation')
584 call glb_intrp_find%init_dofs(this%pe_size)
585 call send_pe_find%init()
586 call recv_pe_find%init()
587 call mpi_barrier(this%comm)
588 time_start = mpi_wtime()
589 write(log_buf, '(A)') 'Global interpolation, finding points'
590 call neko_log%message(log_buf)
591 ! Find pe candidates that the points i want may be at
592 ! Add number to n_points_pe_local
593 !Working arrays
594 this%n_points_pe = 0
595 call this%pe_finder%find_batch(this%xyz, this%n_points, &
596 this%points_at_pe, this%n_points_pe)
597 call mpi_barrier(this%comm)
598 time1 = mpi_wtime()
599 write(log_buf, '(A,E15.7)') &
600 'Found PE candidates time since start of findpts (s):', &
601 time1 - time_start
602 call neko_log%message(log_buf)
603
604 !Send number of points I want to candidates
605 ! n_points_local -> how many points might be at this rank
606 ! n_points_pe_local -> how many points local on this rank that other pes
607 ! might want
608 this%n_points_pe_local = 0
609 this%n_points_local = 0
610 call mpi_reduce_scatter_block(this%n_points_pe, this%n_points_local, &
611 1, mpi_integer, mpi_sum, this%comm, ierr)
612 call mpi_alltoall(this%n_points_pe, 1, mpi_integer,&
613 this%n_points_pe_local, 1, mpi_integer, this%comm, ierr)
614
615 !Set up offset arrays
616 this%n_points_offset_pe_local(0) = 0
617 this%n_points_offset_pe(0) = 0
618 do i = 1, (this%pe_size - 1)
619 this%n_points_offset_pe_local(i) = this%n_points_pe_local(i-1)&
620 + this%n_points_offset_pe_local(i-1)
621 this%n_points_offset_pe(i) = this%n_points_pe(i-1)&
622 + this%n_points_offset_pe(i-1)
623 end do
624 do i = 0, (this%pe_size-1)
625 if (this%n_points_pe(i) .gt. 0) then
626 call send_pe_find%push(i)
627 point_ids => this%points_at_pe(i)%array()
628 do j = 1, this%n_points_pe(i)
629 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+1)
630 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+2)
631 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j)-1)+3)
632 end do
633 end if
634 if (this%n_points_pe_local(i) .gt. 0) then
635 call recv_pe_find%push(i)
636 do j = 1, this%n_points_pe_local(i)
637 call glb_intrp_find%recv_dof(i)%push(3*(j + &
638 this%n_points_offset_pe_local(i) - 1) + 1)
639 call glb_intrp_find%recv_dof(i)%push(3*(j + &
640 this%n_points_offset_pe_local(i) - 1) + 2)
641 call glb_intrp_find%recv_dof(i)%push(3*(j + &
642 this%n_points_offset_pe_local(i) - 1) + 3)
643 end do
644 end if
645 end do
646
647
648
649 call glb_intrp_find%init(send_pe_find, recv_pe_find, this%comm)
650
651 call glb_intrp_find_back%init_dofs(this%pe_size)
652 ii = 0
653 do i = 0, (this%pe_size-1)
654 send_recv => glb_intrp_find%recv_dof(i)%array()
655 do j = 1, glb_intrp_find%recv_dof(i)%size()
656 call glb_intrp_find_back%send_dof(i)%push(send_recv(j))
657 end do
658 send_recv => glb_intrp_find%send_dof(i)%array()
659 do j = 1, glb_intrp_find%send_dof(i)%size()
660 ii = ii + 1
661 call glb_intrp_find_back%recv_dof(i)%push(ii)
662 end do
663 end do
664
665 call glb_intrp_find_back%init(recv_pe_find, send_pe_find, this%comm)
666
667
668 if (allocated(this%xyz_local)) then
669 deallocate(this%xyz_local)
670 end if
671 allocate(this%xyz_local(3, this%n_points_local))
672 call glb_intrp_find%sendrecv(this%xyz, this%xyz_local, this%n_points*3, &
673 this%n_points_local*3)
674
675 call mpi_barrier(this%comm)
676 time1 = mpi_wtime()
677 write(log_buf, '(A,E15.7)') &
678 'Sent to points to PE candidates, time since start of ' &
679 // 'find_points (s):', time1 - time_start
680 call neko_log%message(log_buf)
681
682 !Okay, now we need to find the rst...
683 call all_el_candidates%init()
684
685 if (allocated(n_el_cands)) then
686 deallocate(n_el_cands)
687 end if
688
689 allocate(n_el_cands(this%n_points_local))
691 call this%el_finder%find_batch(this%xyz_local, this%n_points_local, &
692 all_el_candidates, n_el_cands)
693
694 n_point_cand = all_el_candidates%size()
695 if (n_point_cand .gt. 1e8) then
696 print *,'Warning, many point candidates on rank', this%pe_rank, &
697 'cands:', n_point_cand, &
698 'Consider increasing number of ranks'
699 end if
700 call x_t%init(n_point_cand)
701 call y_t%init(n_point_cand)
702 call z_t%init(n_point_cand)
703 ii = 0
705 do i = 1 , this%n_points_local
706 do j = 1, n_el_cands(i)
707 ii = ii + 1
708 x_t%x(ii) = this%xyz_local(1,i)
709 y_t%x(ii) = this%xyz_local(2,i)
710 z_t%x(ii) = this%xyz_local(3,i)
711 end do
712 end do
713
714 call mpi_barrier(this%comm)
715 time1 = mpi_wtime()
716 write(log_buf, '(A,E15.7)') &
717 'Element candidates found, now time for finding rst,time ' // &
718 'since start of find_points (s):', time1 - time_start
719 call neko_log%message(log_buf)
720 call rst_local_cand%init(3,n_point_cand)
721 call resx%init(n_point_cand)
722 call resy%init(n_point_cand)
723 call resz%init(n_point_cand)
724
725 ! Find rst within all element candidates for target xyz (x_t, y_t, z_t)
726 call mpi_barrier(this%comm)
727 time1 = mpi_wtime()
728 el_cands => all_el_candidates%array()
729 if (neko_bcknd_device .eq. 1) then
730 call x_t%copy_from(host_to_device,.false.)
731
732 call y_t%copy_from(host_to_device,.false.)
733
734 call z_t%copy_from(host_to_device,.false.)
735 call device_map(el_cands, el_cands_d,n_point_cand)
736 call device_memcpy(el_cands, el_cands_d,n_point_cand, &
737 host_to_device, .true.)
738 end if
739
740 call this%rst_finder%find(rst_local_cand, &
741 x_t, y_t, z_t, &
742 el_cands, n_point_cand, &
743 resx, resy, resz)
744 if (neko_bcknd_device .eq. 1) then
745 call rst_local_cand%copy_from(device_to_host,.false.)
746 call resx%copy_from(device_to_host,.false.)
747 call resy%copy_from(device_to_host,.false.)
748 call resz%copy_from(device_to_host,.true.)
749 call device_deassociate(el_cands)
750 call device_free(el_cands_d)
751 end if
752 call mpi_barrier(this%comm)
753
754 time2 = mpi_wtime()
755 write(log_buf, '(A,E15.7)') &
756 'Found rst with Newton iteration, time (s):', time2-time1
757 call neko_log%message(log_buf)
758
759 write(log_buf, '(A)') &
760 'Checking validity of points and choosing best candidates.'
761 call neko_log%message(log_buf)
762 call mpi_barrier(this%comm,ierr)
763
764 if (allocated(this%rst_local)) deallocate(this%rst_local)
765 if (allocated(this%el_owner0_local)) deallocate(this%el_owner0_local)
766 allocate(this%rst_local(3,this%n_points_local))
767 allocate(this%el_owner0_local(this%n_points_local))
768 ! Choose the best candidate at this rank
769 ii = 0
770 do i = 1 , this%n_points_local
771 this%xyz_local(1,i) = 10.0
772 this%xyz_local(2,i) = 10.0
773 this%xyz_local(3,i) = 10.0
774 this%rst_local(1,i) = 10.0
775 this%rst_local(2,i) = 10.0
776 this%rst_local(3,i) = 10.0
777 this%el_owner0_local(i) = -1
778 do j = 1, n_el_cands(i)
779 ii = ii + 1
780 if (rst_cmp(this%rst_local(:,i), rst_local_cand%x(:,ii),&
781 this%xyz_local(:,i), [resx%x(ii),resy%x(ii),resz%x(ii)], &
782 this%padding)) then
783 this%rst_local(1,i) = rst_local_cand%x(1,ii)
784 this%rst_local(2,i) = rst_local_cand%x(2,ii)
785
786 this%rst_local(3,i) = rst_local_cand%x(3,ii)
787 this%xyz_local(1,i) = resx%x(ii)
788 this%xyz_local(2,i) = resy%x(ii)
789 this%xyz_local(3,i) = resz%x(ii)
790 this%el_owner0_local(i) = el_cands(ii)
791 end if
792 ! if (this%pe_rank .eq. 0) print *,i, this%rst_local(:,i), &
793 ! this%xyz_local(:,i), this%el_owner0_local(i)
794 end do
795 end do
796 call res%init(3,this%n_points)
797 n_glb_point_cand = sum(this%n_points_pe)
798 if (allocated(rst_results)) deallocate(rst_results)
799 if (allocated(res_results)) deallocate(res_results)
800 if (allocated(el_owner_results)) deallocate(el_owner_results)
801 allocate(rst_results(3,n_glb_point_cand))
802 allocate(res_results(3,n_glb_point_cand))
803 allocate(el_owner_results(n_glb_point_cand))
804 res = 1e2_rp
805 this%rst = 1e2
806 this%pe_owner = -1
807 this%el_owner0 = -1
808 call glb_intrp_find_back%sendrecv(this%xyz_local, res_results, &
809 this%n_points_local*3, n_glb_point_cand*3)
810 call glb_intrp_find_back%sendrecv(this%rst_local, rst_results, &
811 this%n_points_local*3, n_glb_point_cand*3)
812 do i = 1, size(glb_intrp_find_back%send_pe)
813 rank = glb_intrp_find_back%send_pe(i)
814 call mpi_isend(this%el_owner0_local( &
815 this%n_points_offset_pe_local(rank) + 1), &
816 this%n_points_pe_local(rank), &
817 mpi_integer, rank, 0, &
818 this%comm, glb_intrp_find_back%send_buf(i)%request, ierr)
819 glb_intrp_find_back%send_buf(i)%flag = .false.
820 end do
821 do i = 1, size(glb_intrp_find_back%recv_pe)
822 rank = glb_intrp_find_back%recv_pe(i)
823 call mpi_irecv(el_owner_results(this%n_points_offset_pe(rank)+1),&
824 this%n_points_pe(rank), &
825 mpi_integer, rank, 0, &
826 this%comm, glb_intrp_find_back%recv_buf(i)%request, ierr)
827 glb_intrp_find_back%recv_buf(i)%flag = .false.
828 end do
829 call glb_intrp_find_back%nbwait_no_op()
830 ii = 0
831 do i = 1, size(glb_intrp_find_back%recv_pe)
832 point_ids => this%points_at_pe(glb_intrp_find_back%recv_pe(i))%array()
833 do j = 1, this%n_points_pe(glb_intrp_find_back%recv_pe(i))
834 point_id = point_ids(j)
835 ii = ii + 1
836 if (rst_cmp(this%rst(:,point_id), rst_results(:,ii), &
837 res%x(:,point_id), res_results(:,ii), this%padding) .or. &
838 this%pe_owner(point_ids(j)) .eq. -1 ) then
839 this%rst(:,point_ids(j)) = rst_results(:,ii)
840 res%x(:,point_ids(j)) = res_results(:,ii)
841 this%pe_owner(point_ids(j)) = glb_intrp_find_back%recv_pe(i)
842 this%el_owner0(point_ids(j)) = el_owner_results(ii)
843 end if
844 ! if (this%pe_rank .eq. 0) print *,point_id, &
845 !this%rst(:,point_ids(j)),res%x(:,point_ids(j)),
846 ! this%el_owner0(point_ids(j))
847 end do
848 end do
849
850 !OK, now I know the correct rst values of the points I want We now
851 !send the correct rsts to the correct rank (so a point only
852 !belongs to one rank)
853 do i = 0, this%pe_size-1
854 call this%points_at_pe(i)%clear()
855 this%n_points_pe(i) = 0
856 end do
857
858 do i = 1, this%n_points
859 stupid_intent = i
860 if (this%pe_owner(i) .eq. -1 .or. this%el_owner0(i) .eq. -1) then
861 print *, 'No owning rank found for',&
862 ' point ', stupid_intent, ' with coords', this%xyz(:,i), &
863 ' Interpolation will always yield 0.0. Try increase padding.'
864 else
865 call this%points_at_pe(this%pe_owner(i))%push(stupid_intent)
866 this%n_points_pe(this%pe_owner(i)) = &
867 this%n_points_pe(this%pe_owner(i)) + 1
868 end if
869 end do
870 call mpi_reduce_scatter_block(this%n_points_pe, this%n_points_local, 1, &
871 mpi_integer, mpi_sum, this%comm, ierr)
872 call mpi_alltoall(this%n_points_pe, 1, mpi_integer, &
873 this%n_points_pe_local, 1, mpi_integer, this%comm, ierr)
874 this%n_points_offset_pe_local(0) = 0
875 this%n_points_offset_pe(0) = 0
876 do i = 1, (this%pe_size - 1)
877 this%n_points_offset_pe_local(i) = this%n_points_pe_local(i-1)&
878 + this%n_points_offset_pe_local(i-1)
879 this%n_points_offset_pe(i) = this%n_points_pe(i-1)&
880 + this%n_points_offset_pe(i-1)
881 end do
882 call send_pe_find%free()
883 call recv_pe_find%free()
884 call glb_intrp_find%free()
885 call send_pe_find%init()
886 call recv_pe_find%init()
887 call glb_intrp_find%init_dofs(this%pe_size)
888 !setup comm to send xyz and rst to chosen ranks
889 do i = 0, (this%pe_size-1)
890 if (this%n_points_pe(i) .gt. 0) then
891 call send_pe_find%push(i)
892 point_ids => this%points_at_pe(i)%array()
893 do j = 1, this%n_points_pe(i)
894 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 1)
895 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 2)
896 call glb_intrp_find%send_dof(i)%push(3*(point_ids(j) - 1) + 3)
897 end do
898 end if
899 if (this%n_points_pe_local(i) .gt. 0) then
900 call recv_pe_find%push(i)
901 do j = 1, this%n_points_pe_local(i)
902 call glb_intrp_find%recv_dof(i)%push(3*(j + &
903 this%n_points_offset_pe_local(i) - 1) + 1)
904 call glb_intrp_find%recv_dof(i)%push(3*(j + &
905 this%n_points_offset_pe_local(i) - 1) + 2)
906 call glb_intrp_find%recv_dof(i)%push(3*(j + &
907 this%n_points_offset_pe_local(i) - 1) + 3)
908 end do
909 end if
910 end do
911
912
913 call glb_intrp_find%init(send_pe_find, recv_pe_find, this%comm)
914 call glb_intrp_find%sendrecv(this%xyz, this%xyz_local, this%n_points*3, &
915 this%n_points_local*3)
916 call glb_intrp_find%sendrecv(this%rst, this%rst_local, this%n_points*3, &
917 this%n_points_local*3)
918 ii = 0
919 do i = 1, size(glb_intrp_find%send_pe)
920 rank = glb_intrp_find%send_pe(i)
921 point_ids => this%points_at_pe(rank)%array()
922 do j = 1, this%n_points_pe(rank)
923 ii = ii + 1
924 el_owner_results(ii) = this%el_owner0(point_ids(j))
925 end do
926 call mpi_isend(el_owner_results(this%n_points_offset_pe(rank) + 1),&
927 this%n_points_pe(rank), &
928 mpi_integer, rank, 0, &
929 this%comm, glb_intrp_find%send_buf(i)%request, ierr)
930 glb_intrp_find%send_buf(i)%flag = .false.
931 end do
932 do i = 1, size(glb_intrp_find%recv_pe)
933 rank = glb_intrp_find%recv_pe(i)
934 call mpi_irecv(this%el_owner0_local( &
935 this%n_points_offset_pe_local(rank) + 1), &
936 this%n_points_pe_local(rank), &
937 mpi_integer, rank, 0, &
938 this%comm, glb_intrp_find%recv_buf(i)%request, ierr)
939 glb_intrp_find%recv_buf(i)%flag = .false.
940 end do
941 call glb_intrp_find%nbwait_no_op()
942
943 call glb_intrp_find%free()
944
945 !Set up final way of doing communication
946 call send_pe%init()
947 call recv_pe%init()
948 call this%glb_intrp_comm%init_dofs(this%pe_size)
949 do i = 0, (this%pe_size-1)
950 if (this%n_points_pe(i) .gt. 0) then
951 call recv_pe%push(i)
952 point_ids => this%points_at_pe(i)%array()
953 do j = 1, this%n_points_pe(i)
954 call this%glb_intrp_comm%recv_dof(i)%push(point_ids(j))
955 end do
956 end if
957 if (this%n_points_pe_local(i) .gt. 0) then
958 call send_pe%push(i)
959 do j = 1, this%n_points_pe_local(i)
960 call this%glb_intrp_comm%send_dof(i)%push(j + &
961 this%n_points_offset_pe_local(i))
962 end do
963 end if
964 end do
965 call this%glb_intrp_comm%init(send_pe, recv_pe,this%comm)
966
967 !Initialize working arrays for evaluation
968 call this%temp_local%init(this%n_points_local)
969 call this%temp%init(this%n_points)
970
971 !Initialize interpolator for local interpolation
972 call this%local_interp%init(this%Xh, this%rst_local, &
973 this%n_points_local)
974
975
976 if (neko_bcknd_device .eq. 1) then
977 call device_memcpy(this%el_owner0, this%el_owner0_d, &
978 this%n_points, host_to_device, sync = .true.)
979 call device_map(this%el_owner0_local, this%el_owner0_local_d, &
980 this%n_points_local)
981 call device_memcpy(this%el_owner0_local, this%el_owner0_local_d, &
982 this%n_points_local, host_to_device, sync = .true.)
983 end if
984
985 call this%check_points(this%x%x,this%y%x, this%z%x)
986
987 !Free stuff
988 call send_pe%free()
989 call recv_pe%free()
990 call glb_intrp_find_back%free()
991 call send_pe_find%free()
992 call recv_pe_find%free()
993 call x_t%free()
994 call y_t%free()
995 call z_t%free()
996 call rst_local_cand%free()
997 call resx%free()
998 call resy%free()
999 call resz%free()
1000 call res%free()
1001 call all_el_candidates%free()
1002
1003 if (allocated(n_el_cands)) deallocate(n_el_cands)
1004 if (allocated(rst_results)) deallocate(rst_results)
1005 if (allocated(res_results)) deallocate(res_results)
1006 if (allocated(el_owner_results)) deallocate(el_owner_results)
1007 call mpi_barrier(this%comm, ierr)
1008 time2 = mpi_wtime()
1009 write(log_buf, '(A,E15.7)') 'Global interpolation find points ' // &
1010 'done, time (s):', time2-time_start
1011 call neko_log%message(log_buf)
1012 call neko_log%end_section()
1013 call neko_log%newline()
1015
1019 subroutine global_interpolation_check_points(this, x, y, z)
1020 class(global_interpolation_t), target, intent(inout) :: this
1021 real(kind=rp), intent(inout) :: x(:)
1022 real(kind=rp), intent(inout) :: y(:)
1023 real(kind=rp), intent(inout) :: z(:)
1024 integer :: i, j
1025 character(len=8000) :: log_buf
1026 real(kind=rp) :: xdiff, ydiff, zdiff
1027 logical :: isdiff
1028 type(vector_t) :: x_check, y_check, z_check
1029
1030 call x_check%init(this%n_points)
1031 call y_check%init(this%n_points)
1032 call z_check%init(this%n_points)
1033 call this%evaluate(x_check%x, x, on_host = .true.)
1034 call this%evaluate(y_check%x, y, on_host = .true.)
1035 call this%evaluate(z_check%x, z, on_host = .true.)
1036 write(log_buf, '(A)') 'Checking validity of points.'
1037 call neko_log%message(log_buf)
1038 j = 0
1039 do i = 1 , this%n_points
1040 ! Check validity of points
1041 isdiff = .false.
1042 xdiff = x_check%x(i)-this%xyz(1,i)
1043 ydiff = y_check%x(i)-this%xyz(2,i)
1044 zdiff = z_check%x(i)-this%xyz(3,i)
1045 isdiff = norm2(real([xdiff,ydiff,zdiff],xp)) > this%tolerance
1046 if (isdiff) then
1047 write(*,*) 'Point ', i,'at rank ', this%pe_rank, &
1048 'with coordinates: ', &
1049 this%xyz(1, i), this%xyz(2, i), this%xyz(3, i), &
1050 'Differ from interpolated coords: ', &
1051 x_check%x(i), y_check%x(i), z_check%x(i), &
1052 'Actual difference: ', &
1053 xdiff, ydiff, zdiff, norm2(real([xdiff,ydiff,zdiff],xp)),&
1054 'Process, element: ', &
1055 this%pe_owner(i), this%el_owner0(i)+1, &
1056 'Calculated rst: ', &
1057 this%rst(1,i), this%rst(2,i), this%rst(3,i)
1058 j = j + 1
1059 end if
1060 end do
1061 call x_check%free()
1062 call y_check%free()
1063 call z_check%free()
1064
1066
1074 subroutine global_interpolation_find_coords(this, x, y, z, n_points)
1075 class(global_interpolation_t), intent(inout) :: this
1076 integer :: n_points
1077 !!Perhaps this should be kind dp
1078 !!this is to get around that x,y,z often is 4 dimensional...
1079 !!Should maybe add interface for 1d aswell
1080 real(kind=rp) :: x(n_points,1,1,1)
1081 real(kind=rp) :: y(n_points,1,1,1)
1082 real(kind=rp) :: z(n_points,1,1,1)
1083 integer :: i
1084
1085 call this%free_points()
1086 call this%free_points_local()
1087
1088 this%n_points = n_points
1089
1091
1092 !Deepcopy of coordinates
1093 do i = 1, n_points
1094 this%xyz(1, i) = x(i,1,1,1)
1095 this%xyz(2, i) = y(i,1,1,1)
1096 this%xyz(3, i) = z(i,1,1,1)
1097 end do
1098
1100
1109 subroutine global_interpolation_find_coords1d(this, x, y, z, n_points)
1110 class(global_interpolation_t), intent(inout) :: this
1111 integer :: n_points
1112 real(kind=rp) :: x(n_points)
1113 real(kind=rp) :: y(n_points)
1114 real(kind=rp) :: z(n_points)
1115 integer :: i
1116
1117 call this%free_points()
1118 call this%free_points_local()
1119
1120 this%n_points = n_points
1121
1123
1124 !Deepcopy of coordinates
1125 do i = 1, n_points
1126 this%xyz(1, i) = x(i)
1127 this%xyz(2, i) = y(i)
1128 this%xyz(3, i) = z(i)
1129 end do
1130
1132
1134
1135
1137 class(global_interpolation_t) :: this
1138
1139 if (allocated(this%xyz)) deallocate(this%xyz)
1140 if (allocated(this%rst)) deallocate(this%rst)
1141 if (allocated(this%pe_owner)) deallocate(this%pe_owner)
1142 if (allocated(this%el_owner0)) deallocate(this%el_owner0)
1143
1144 allocate(this%pe_owner(this%n_points))
1145 allocate(this%el_owner0(this%n_points))
1146 allocate(this%xyz(3, this%n_points))
1147 allocate(this%rst(3, this%n_points))
1148 if (neko_bcknd_device .eq. 1) then
1149 call device_map(this%el_owner0, this%el_owner0_d, this%n_points)
1150 end if
1151
1153
1160 subroutine global_interpolation_find_xyz(this, xyz, n_points)
1161 class(global_interpolation_t), intent(inout) :: this
1162 integer, intent(in) :: n_points
1163 !!Perhaps this should be kind dp
1164 real(kind=rp), intent(inout) :: xyz(3, n_points)
1165
1166
1167 call this%free_points()
1168 call this%free_points_local()
1169
1170 this%n_points = n_points
1171
1173
1175 call copy(this%xyz, xyz, 3 * n_points)
1176
1178
1179 end subroutine global_interpolation_find_xyz
1180
1189 subroutine global_interpolation_find_and_redist(this, xyz, n_points)
1190 class(global_interpolation_t), intent(inout) :: this
1191 integer, intent(inout) :: n_points
1192 !!Perhaps this should be kind dp
1193 real(kind=rp), allocatable, intent(inout) :: xyz(:,:)
1194
1195
1196 call this%free_points()
1197 call this%free_points_local()
1198
1199
1200 this%n_points = n_points
1202
1204 call copy(this%xyz, xyz, 3 * n_points)
1205
1207 call this%free_points()
1208 this%n_points = this%n_points_local
1209 n_points = this%n_points_local
1211 if (allocated(xyz)) then
1212 deallocate(xyz)
1213 end if
1214 allocate(xyz(3,n_points))
1215
1216 call copy(xyz, this%xyz_local, 3*n_points)
1217 call copy(this%rst, this%rst_local, 3*n_points)
1218 call copy(this%xyz, this%xyz_local, 3*n_points)
1219 this%pe_owner = this%pe_rank
1220 this%el_owner0 = this%el_owner0_local
1221 if (neko_bcknd_device .eq. 1) then
1222 call device_memcpy(this%el_owner0, this%el_owner0_d, &
1223 this%n_points, host_to_device, sync = .true.)
1224 end if
1225 this%all_points_local = .true.
1226
1228
1235 subroutine global_interpolation_evaluate_masked(this, interp_values, &
1236 field, mask, on_host)
1237 class(global_interpolation_t), target, intent(inout) :: this
1238 real(kind=rp), intent(inout), target :: interp_values(this%n_points)
1239 real(kind=rp), intent(inout), target :: field(this%n_dof)
1240 type(mask_t), intent(in) :: mask
1241 logical, intent(in) :: on_host
1242
1243 call vector_masked_gather_copy(this%masked_field, field, mask, this%n_dof)
1244 call this%evaluate(interp_values, this%masked_field%x, on_host)
1245
1247
1252 subroutine global_interpolation_evaluate(this, interp_values, field, on_host)
1253 class(global_interpolation_t), target, intent(inout) :: this
1254 real(kind=rp), intent(inout), target :: interp_values(this%n_points)
1255 real(kind=rp), intent(inout), target :: field(this%nelv*this%Xh%lxyz)
1256 logical, intent(in) :: on_host
1257 type(c_ptr) :: interp_d
1258
1259 if (.not. this%all_points_local) then
1260 call this%local_interp%evaluate(this%temp_local%x, &
1261 this%el_owner0_local, field, this%nelv, on_host)
1262 if (neko_bcknd_device .eq. 1 .and. .not. on_host) then
1263 call device_memcpy(this%temp_local%x, this%temp_local%x_d, &
1264 this%n_points_local, device_to_host, .true.)
1265 end if
1266 interp_values = 0.0_rp
1267 call this%glb_intrp_comm%sendrecv(this%temp_local%x, interp_values, &
1268 this%n_points_local, this%n_points)
1269 if (neko_bcknd_device .eq. 1 .and. .not. on_host) then
1270 interp_d = device_get_ptr(interp_values)
1271 call device_memcpy(interp_values, interp_d, &
1272 this%n_points, host_to_device, .false.)
1273 end if
1274 else
1275 call this%local_interp%evaluate(interp_values, this%el_owner0_local, &
1276 field, this%nelv, on_host)
1277 end if
1278
1279 end subroutine global_interpolation_evaluate
1280
1281
1292 function rst_cmp(rst1, rst2,res1, res2, tol) result(rst2_better)
1293 real(kind=rp) :: rst1(3), res1(3)
1294 real(kind=rp) :: rst2(3), res2(3)
1295 real(kind=dp) :: tol
1296 logical :: rst2_better
1297 !If rst1 is invalid and rst2 is valid, take rst2
1298 ! If both invalidl, take smallest residual
1299 if (abs(rst1(1)) .gt. 1.0_xp+tol .or. &
1300 abs(rst1(2)) .gt. 1.0_xp+tol .or. &
1301 abs(rst1(3)) .gt. 1.0_xp+tol) then
1302 if (abs(rst2(1)) .le. 1.0_xp+tol .and. &
1303 abs(rst2(2)) .le. 1.0_xp+tol .and. &
1304 abs(rst2(3)) .le. 1.0_xp+tol) then
1305 rst2_better = .true.
1306 else
1307 rst2_better = (norm2(real(res2,xp)) .lt. norm2(real(res1,xp)))
1308 end if
1309 else
1311 rst2_better = (norm2(real(res2,xp)) .lt. norm2(real(res1,xp)) .and.&
1312 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)
1315 end if
1316 end function rst_cmp
1317
1318end 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:250
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: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 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