Neko  0.8.99
A portable framework for high-order spectral element flow simulations
global_interpolation.F90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2023, 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 !
39  use num_types, only: rp
40  use space, only: space_t
41  use dofmap, only: dofmap_t
42  use mesh, only: mesh_t
43  use logger, only: neko_log, log_size
44  use utils, only: neko_error, neko_warning
46  use comm
47  use math, only: copy
48  use neko_mpi_types
49  use, intrinsic :: iso_c_binding
50  implicit none
51  private
53  type, public :: global_interpolation_t
55  type(dofmap_t), pointer :: dof
57  type(mesh_t), pointer :: mesh
59  type(space_t), pointer :: xh
61  type(local_interpolator_t) :: local_interp
63  logical :: all_points_local = .false.
64  !! Gslib handle
65  !! @note: Remove when we remove gslib
66  integer :: gs_handle
67  logical :: gs_init = .false.
70  integer :: n_points
73  real(kind=rp), allocatable :: xyz(:,:)
75  integer, allocatable :: proc_owner(:)
77  integer, allocatable :: el_owner(:)
78  type(c_ptr) :: el_owner_d = c_null_ptr
81  real(kind=rp), allocatable :: rst(:,:)
84  real(kind=rp), allocatable :: dist2(:)
86  integer, allocatable :: error_code(:)
88  real(kind=rp) :: tol = 5e-13
89  contains
91  procedure, pass(this) :: init => global_interpolation_init
93  procedure, pass(this) :: free => global_interpolation_free
95  procedure, pass(this) :: free_points => global_interpolation_free_points
96  procedure, pass(this) :: find_points_and_redist => &
101  procedure, pass(this) :: find_points_coords => global_interpolation_find_coords
102  procedure, pass(this) :: find_points_xyz => global_interpolation_find_xyz
103  generic :: find_points => find_points_xyz, find_points_coords
105  procedure, pass(this) :: evaluate => global_interpolation_evaluate
107  end type global_interpolation_t
108 
109 
110 contains
114  subroutine global_interpolation_init(this, dof, tol)
115  class(global_interpolation_t), intent(inout) :: this
116  type(dofmap_t), target :: dof
117  real(kind=rp), optional :: tol
118  integer :: lx, ly, lz, nelv, max_pts_per_iter
119 
120  this%dof => dof
121  this%Xh => dof%Xh
122  this%mesh => dof%msh
123  if(present(tol)) this%tol = tol
124 
125 #ifdef HAVE_GSLIB
126 
127  lx = this%Xh%lx
128  ly = this%Xh%ly
129  lz = this%Xh%lz
130  nelv = this%mesh%nelv
131  !Number of points to iterate on simultaneosuly
132  max_pts_per_iter = 128
133 
134  call fgslib_findpts_setup(this%gs_handle, &
135  neko_comm, pe_size, &
136  this%mesh%gdim, &
137  dof%x, dof%y, dof%z, & ! Physical nodal values
138  lx, ly, lz, nelv, & ! Mesh dimensions
139  2*lx, 2*ly, 2*lz, & ! Mesh size for bounding box computation
140  0.01, & ! relative size to expand bounding boxes by
141  lx*ly*lz*nelv, lx*ly*lz*nelv, & ! local/global hash mesh sizes
142  max_pts_per_iter, this%tol)
143  this%gs_init = .true.
144 #else
145  call neko_error('Neko needs to be built with GSLIB support')
146 #endif
147 
148  end subroutine global_interpolation_init
149 
150 
152  subroutine global_interpolation_free(this)
153  class(global_interpolation_t), intent(inout) :: this
154 
155  nullify(this%mesh)
156  nullify(this%dof)
157  nullify(this%Xh)
158 
159  call this%free_points()
160 
161 #ifdef HAVE_GSLIB
162  if (this%gs_init) then
163  call fgslib_findpts_free(this%gs_handle)
164  this%gs_init = .false.
165  end if
166 #endif
167 
168  end subroutine global_interpolation_free
169 
172  class(global_interpolation_t), intent(inout) :: this
173 
174  this%n_points = 0
175  this%all_points_local = .false.
176 
177  if (allocated(this%xyz)) deallocate(this%xyz)
178  if (allocated(this%rst)) deallocate(this%rst)
179  if (allocated(this%proc_owner)) deallocate(this%proc_owner)
180  if (allocated(this%el_owner)) deallocate(this%el_owner)
181  if (allocated(this%dist2)) deallocate(this%dist2)
182  if (allocated(this%error_code)) deallocate(this%error_code)
183 
184  if (c_associated(this%el_owner_d)) then
185  call device_free(this%el_owner_d)
186  end if
187 
188  end subroutine global_interpolation_free_points
189 
192  class(global_interpolation_t), intent(inout) :: this
193  !!Perhaps this should be kind dp
194  real(kind=rp) :: xdiff, ydiff, zdiff
195  character(len=8000) :: log_buf
196  real(kind=rp), allocatable :: x_check(:)
197  real(kind=rp), allocatable :: y_check(:)
198  real(kind=rp), allocatable :: z_check(:)
199  logical :: isdiff
200  integer :: i
201 
202 
203 #ifdef HAVE_GSLIB
204 
205  ! gslib find points, which element they belong, to process etc.
206  call fgslib_findpts(this%gs_handle, &
207  this%error_code, 1, &
208  this%proc_owner, 1, &
209  this%el_owner, 1, &
210  this%rst, this%mesh%gdim, &
211  this%dist2, 1, &
212  this%xyz(1,1), this%mesh%gdim, &
213  this%xyz(2,1), this%mesh%gdim, &
214  this%xyz(3,1), this%mesh%gdim, this%n_points)
215 
216  do i=1,this%n_points
217 
218  !
219  ! Check validity of points
220  !
221  if (this%error_code(i) .eq. 1) then
222  if (this%dist2(i) .gt. this%tol) then
223  write(*,*) 'Point with coords: ',&
224  this%xyz(1,i),&
225  this%xyz(2,i),&
226  this%xyz(3,i),&
227  'Did not converge to tol. Absolute differences squared: ',&
228  this%dist2(i), 'PE rank', pe_rank
229  end if
230  end if
231 
232  if (this%error_code(i) .eq. 2) &
233  write(*,*) 'Point with coords: ',&
234  this%xyz(1,i), this%xyz(2,i), this%xyz(3,i),&
235  'Outside the mesh!',&
236  ' Interpolation on these points will return 0.0. dist2: ', &
237  this%dist2(i),&
238  'el_owner, rst coords, pe: ',&
239  this%el_owner(i), this%rst(1,i), this%rst(2,i), this%rst(3,i), pe_rank
240 
241  end do
242 
243  allocate(x_check(this%n_points))
244  allocate(y_check(this%n_points))
245  allocate(z_check(this%n_points))
246 
247  call fgslib_findpts_eval(this%gs_handle, x_check, &
248  1, this%error_code, 1, &
249  this%proc_owner, 1, this%el_owner, 1, &
250  this%rst, this%mesh%gdim, &
251  this%n_points, this%dof%x)
252 
253  call fgslib_findpts_eval(this%gs_handle, y_check, &
254  1, this%error_code, 1, &
255  this%proc_owner, 1, this%el_owner, 1, &
256  this%rst, this%mesh%gdim, &
257  this%n_points, this%dof%y)
258 
259  call fgslib_findpts_eval(this%gs_handle, z_check, &
260  1, this%error_code, 1, &
261  this%proc_owner, 1, this%el_owner, 1, &
262  this%rst, this%mesh%gdim, &
263  this%n_points, this%dof%z)
264 
265 
266  do i=1,this%n_points
267 
268  !
269  ! Check validity of points
270  !
271  isdiff = .false.
272  xdiff = (x_check(i)-this%xyz(1,i))**2
273  if ( xdiff .gt. this%tol) isdiff = .true.
274  ydiff = (y_check(i)-this%xyz(2,i))**2
275  if ( ydiff .gt. this%tol) isdiff = .true.
276  zdiff = (z_check(i)-this%xyz(3,i))**2
277  if ( zdiff .gt. this%tol) isdiff = .true.
278 
279  if (isdiff) then
280  write(*,*) 'Points with coords: ',&
281  this%xyz(1,i),this%xyz(2,i),this%xyz(3,i), &
282  'Differ from interpolated coords: ',&
283  x_check(i), y_check(i), z_check(i),&
284  'Distance squared: ',&
285  xdiff, ydiff, zdiff
286  end if
287 
288  end do
289 
290  deallocate(x_check)
291  deallocate(y_check)
292  deallocate(z_check)
293 
294  if (neko_bcknd_device .eq. 1) then
295  call device_memcpy(this%el_owner, this%el_owner_d, &
296  this%n_points, host_to_device, sync = .true.)
297  end if
298 #else
299  call neko_error('Neko needs to be built with GSLIB support')
300 #endif
301  end subroutine global_interpolation_find_common
302 
315  subroutine global_interpolation_find_coords(this, x, y, z, n_points)
316  class(global_interpolation_t), intent(inout) :: this
317  integer :: n_points
318  !!Perhaps this should be kind dp
319  !!this is to get around that x,y,z often is 4 dimensional...
320  !!Should maybe add interface for 1d aswell
321  real(kind=rp) :: x(n_points,1,1,1)
322  real(kind=rp) :: y(n_points,1,1,1)
323  real(kind=rp) :: z(n_points,1,1,1)
324  integer :: i
325 
326  call this%free_points()
327 
328  this%n_points = n_points
329 
331 
332  do i = 1, n_points
333  this%xyz(1,i) = x(i,1,1,1)
334  this%xyz(2,i) = y(i,1,1,1)
335  this%xyz(3,i) = z(i,1,1,1)
336  end do
337 
339 
340  end subroutine global_interpolation_find_coords
341 
343  class(global_interpolation_t) :: this
344 
345  allocate(this%xyz(3,this%n_points))
346  allocate(this%rst(3,this%n_points))
347  allocate(this%proc_owner(this%n_points))
348  allocate(this%el_owner(this%n_points))
349  allocate(this%dist2(this%n_points))
350  allocate(this%error_code(this%n_points))
351 
352  if (neko_bcknd_device .eq. 1) &
353  call device_map(this%el_owner, this%el_owner_d,this%n_points)
354 
356 
367  subroutine global_interpolation_find_xyz(this, xyz, n_points)
368  class(global_interpolation_t), intent(inout) :: this
369  integer, intent(in) :: n_points
370  !!Perhaps this should be kind dp
371  real(kind=rp), intent(inout) :: xyz(3,n_points)
372 
373 
374  call this%free_points()
375 
376  this%n_points = n_points
377 
379 
381  call copy(this%xyz,xyz,3*n_points)
382 
384 
385  end subroutine global_interpolation_find_xyz
386 
397  subroutine global_interpolation_find_and_redist(this, xyz, n_points)
398  class(global_interpolation_t), intent(inout) :: this
399  integer, intent(inout) :: n_points
400  !!Perhaps this should be kind dp
401  real(kind=rp), allocatable, intent(inout) :: xyz(:,:)
402  integer :: i
403 
404 
405  call this%free_points()
406 
407  this%n_points = n_points
409 
411  call copy(this%xyz,xyz,3*n_points)
412 
415  call global_interpolation_redist(this)
417 
418  do i = 1, this%n_points
419  if (this%proc_owner(i) .ne. pe_rank) then
420  write(*,*) 'Redistribution failed on rank: ', pe_rank,&
421  'for point with coord: ', &
422  this%xyz(1,i),this%xyz(2,i),this%xyz(3,i)
423  exit
424  end if
425  end do
426 
427  n_points = this%n_points
428  deallocate(xyz)
429  allocate(xyz(3,n_points))
430  call copy(xyz,this%xyz,3*n_points)
431 
432  call this%local_interp%init(this%Xh, this%rst(1,:),&
433  this%rst(2,:), this%rst(3,:), n_points)
434  this%all_points_local = .true.
435 
436 
438 
440  class(global_interpolation_t), intent(inout) :: this
441  integer, allocatable :: n_points_per_pe(:)
442  integer, allocatable :: n_points_from_pe(:)
443  integer, allocatable :: n_point_offset_from_pe(:)
444  real(kind=rp), allocatable :: xyz_send_to_pe(:,:)
445  real(kind=rp), allocatable :: new_xyz(:,:)
446  integer :: i, j, k, ierr, n_new_points, max_n_points_to_send
447 
448  n_new_points = 0
449 
450  allocate(n_points_per_pe(0:(pe_size-1)))
451  allocate(n_points_from_pe(0:(pe_size-1)))
452  n_points_per_pe = 0
453  n_points_from_pe = 0
455  do i = 1,this%n_points
456  n_points_per_pe(this%proc_owner(i)) = n_points_per_pe(this%proc_owner(i)) + 1
457  end do
460  do i = 0,(pe_size-1)
461  call mpi_reduce(n_points_per_pe(i),n_new_points,1,mpi_integer,mpi_sum, i, neko_comm, ierr)
462  call mpi_gather(n_points_per_pe(i), 1, mpi_integer,&
463  n_points_from_pe, 1, mpi_integer, i, neko_comm, ierr)
464  end do
465 
466  allocate(n_point_offset_from_pe(0:(pe_size-1)))
467  n_point_offset_from_pe(0) = 0
468  do i = 1,(pe_size-1)
469  n_point_offset_from_pe(i) = n_points_from_pe(i-1)&
470  + n_point_offset_from_pe(i-1)
471  end do
472 
473  allocate(new_xyz(3,n_new_points))
474  max_n_points_to_send = maxval(n_points_per_pe)
475  allocate(xyz_send_to_pe(3,max_n_points_to_send))
476  do i = 0, (pe_size-1)
478  k = 0
479  do j = 1, this%n_points
480  if (this%proc_owner(j) .eq. i) then
481  k = k + 1
482  xyz_send_to_pe(:,k) = this%xyz(:,j)
483  end if
484  end do
485  if (k .ne. n_points_per_pe(i)) then
486  write(*,*) 'PE: ', pe_rank, ' has k= ', k,&
487  'points for PE:', i,' but should have: ',&
488  n_points_per_pe(i)
489  call neko_error('Error in redistribution of points')
490  end if
491  call mpi_gatherv(xyz_send_to_pe,3*n_points_per_pe(i),&
492  mpi_double_precision, new_xyz,3*n_points_from_pe,&
493  3*n_point_offset_from_pe,&
494  mpi_double_precision, i, neko_comm, ierr)
495 
496  end do
497 
498  call this%free_points()
499 
500  this%n_points = n_new_points
502  call copy(this%xyz,new_xyz,3*n_new_points)
503 
504  deallocate(n_point_offset_from_pe)
505  deallocate(n_points_from_pe)
506  deallocate(n_points_per_pe)
507  deallocate(xyz_send_to_pe)
508 
509  end subroutine global_interpolation_redist
510 
511 
512 
516  subroutine global_interpolation_evaluate(this, interp_values, field)
517  class(global_interpolation_t), intent(inout) :: this
518  real(kind=rp), intent(inout) :: interp_values(this%n_points)
519  real(kind=rp), intent(inout) :: field(this%dof%size())
520 
521 #ifdef HAVE_GSLIB
522  if (.not. this%all_points_local) then
523  call fgslib_findpts_eval(this%gs_handle, interp_values, &
524  1, this%error_code, 1, &
525  this%proc_owner, 1, this%el_owner, 1, &
526  this%rst, this%mesh%gdim, &
527  this%n_points, field)
528  else
529  if (this%n_points .gt. 0) &
530  call this%local_interp%evaluate(interp_values, this%el_owner,&
531  field, this%mesh%nelv)
532  end if
533 #else
534  call neko_error('Neko needs to be built with GSLIB support')
535 #endif
536 
537  end subroutine global_interpolation_evaluate
538 
539 end module global_interpolation
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:26
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
integer pe_size
MPI size of communicator.
Definition: comm.F90:29
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Defines a field.
Definition: field.f90:34
Implements global_interpolation given a dofmap.
subroutine global_interpolation_free(this)
Destructor.
subroutine global_interpolation_free_points(this)
Destructor for point arrays.
subroutine global_interpolation_init_point_arrays(this)
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...
subroutine global_interpolation_evaluate(this, interp_values, field)
Evalute the interpolated value in the points given a field on the dofmap.
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_redist(this)
subroutine global_interpolation_init(this, dof, tol)
Evaluate only local points.
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:61
integer, parameter, public log_size
Definition: log.f90:40
Definition: math.f90:60
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:228
Defines a mesh.
Definition: mesh.f90:34
MPI derived types.
Definition: mpi_types.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a function space.
Definition: space.f90:34
Utilities.
Definition: utils.f90:35
subroutine, public neko_warning(warning_msg)
Definition: utils.f90:198
Implements global interpolation for arbitrary points in the domain.
Interpolation on a set of points with known rst coordinates in elements local to this process....
The function space for the SEM solution fields.
Definition: space.f90:62