Neko  0.9.0
A portable framework for high-order spectral element flow simulations
point_interpolator.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-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 !
36  use tensor, only: triple_tensor_product
37  use space, only: space_t, gl, gll
38  use num_types, only: rp
39  use point, only: point_t
40  use math, only: abscmp
41  use fast3d, only: fd_weights_full
42  use utils, only: neko_error
43  use device
44  use device_math, only: device_rzero
46  implicit none
47  private
48 
55  type, public :: point_interpolator_t
57  type(space_t), pointer :: xh => null()
58  contains
60  procedure, pass(this) :: init => point_interpolator_init
62  procedure, pass(this) :: free => point_interpolator_free
64  procedure, pass(this) :: compute_weights => point_interpolator_compute_weights
66  procedure, pass(this) :: point_interpolator_interpolate_scalar
68  procedure, pass(this) :: point_interpolator_interpolate_vector
70  procedure, pass(this) :: point_interpolator_interpolate_jacobian
74  generic :: interpolate => point_interpolator_interpolate_scalar, &
77  generic :: jacobian => point_interpolator_interpolate_jacobian, &
79 
80  end type point_interpolator_t
81 
82 contains
83 
86  subroutine point_interpolator_init(this, Xh)
87  class(point_interpolator_t), intent(inout), target :: this
88  type(space_t), intent(in), target :: Xh
89 
90  if ((xh%t .eq. gl) .or. (xh%t .eq. gll)) then
91  else
92  call neko_error('Unsupported interpolation')
93  end if
94 
95  this%Xh => xh
96 
97  end subroutine point_interpolator_init
98 
100  subroutine point_interpolator_free(this)
101  class(point_interpolator_t), intent(inout) :: this
102 
103  if (associated(this%Xh)) this%Xh => null()
104 
105  end subroutine point_interpolator_free
106 
117  subroutine point_interpolator_compute_weights(this, r, s, t, wr, ws, wt)
118  class(point_interpolator_t), intent(inout) :: this
119  real(kind=rp), intent(in) :: r(:), s(:), t(:)
120  real(kind=rp), intent(inout) :: wr(:,:), ws(:,:), wt(:,:)
121 
122  integer :: N, i, lx
123  lx = this%Xh%lx
124  n = size(r)
125 
126  do i = 1, n
127  call fd_weights_full(r(i), this%Xh%zg(:,1), lx-1, 0, wr(:,i))
128  call fd_weights_full(s(i), this%Xh%zg(:,2), lx-1, 0, ws(:,i))
129  call fd_weights_full(t(i), this%Xh%zg(:,3), lx-1, 0, wt(:,i))
130  end do
131 
133 
139  function point_interpolator_interpolate_scalar(this, rst, X) result(res)
140  class(point_interpolator_t), intent(in) :: this
141  type(point_t), intent(in) :: rst(:)
142  real(kind=rp), intent(inout) :: x(this%Xh%lx, this%Xh%ly, this%Xh%lz)
143  real(kind=rp), allocatable :: res(:)
144 
145  real(kind=rp) :: hr(this%Xh%lx), hs(this%Xh%ly), ht(this%Xh%lz)
146  integer :: lx,ly,lz, i
147  integer :: n
148  lx = this%Xh%lx
149  ly = this%Xh%ly
150  lz = this%Xh%lz
151 
152  n = size(rst)
153  allocate(res(n))
154 
155  !
156  ! Compute weights and perform interpolation for the first point
157  !
158  call fd_weights_full(real(rst(1)%x(1), rp), this%Xh%zg(:,1), lx-1, 0, hr)
159  call fd_weights_full(real(rst(1)%x(2), rp), this%Xh%zg(:,2), ly-1, 0, hs)
160  call fd_weights_full(real(rst(1)%x(3), rp), this%Xh%zg(:,3), lz-1, 0, ht)
161 
162  ! And interpolate!
163  call triple_tensor_product(res(1),x,lx,hr,hs,ht)
164 
165  if (n .eq. 1) return
166 
167  !
168  ! Loop through the rest of the points
169  !
170  do i = 2, n
171 
172  ! If the coordinates are different, then recompute weights
173  if ( .not. abscmp(rst(i)%x(1), rst(i-1)%x(1)) ) then
174  call fd_weights_full(real(rst(i)%x(1), rp), &
175  this%Xh%zg(:,1), lx-1, 0, hr)
176  end if
177  if ( .not. abscmp(rst(i)%x(2), rst(i-1)%x(2)) ) then
178  call fd_weights_full(real(rst(i)%x(2), rp), &
179  this%Xh%zg(:,2), ly-1, 0, hs)
180  end if
181  if ( .not. abscmp(rst(i)%x(3), rst(i-1)%x(3)) ) then
182  call fd_weights_full(real(rst(i)%x(3), rp), &
183  this%Xh%zg(:,3), lz-1, 0, ht)
184  end if
185 
186  ! And interpolate!
187  call triple_tensor_product(res(i), x, lx, hr, hs, ht)
188 
189  end do
190 
192 
200  function point_interpolator_interpolate_vector(this, rst, X, Y, Z) result(res)
201  class(point_interpolator_t), intent(in) :: this
202  type(point_t), intent(in) :: rst(:)
203  real(kind=rp), intent(inout) :: x(this%Xh%lx, this%Xh%ly, this%Xh%lz)
204  real(kind=rp), intent(inout) :: y(this%Xh%lx, this%Xh%ly, this%Xh%lz)
205  real(kind=rp), intent(inout) :: z(this%Xh%lx, this%Xh%ly, this%Xh%lz)
206 
207  type(point_t), allocatable :: res(:)
208  real(kind=rp), allocatable :: tmp(:,:)
209  real(kind=rp) :: hr(this%Xh%lx), hs(this%Xh%ly), ht(this%Xh%lz)
210  integer :: lx,ly,lz, i
211  integer :: n
212  lx = this%Xh%lx
213  ly = this%Xh%ly
214  lz = this%Xh%lz
215 
216  n = size(rst)
217  allocate(res(n))
218  allocate(tmp(3, n))
219 
220  !
221  ! Compute weights and perform interpolation for the first point
222  !
223  call fd_weights_full(real(rst(1)%x(1), rp), this%Xh%zg(:,1), lx-1, 0, hr)
224  call fd_weights_full(real(rst(1)%x(2), rp), this%Xh%zg(:,2), ly-1, 0, hs)
225  call fd_weights_full(real(rst(1)%x(3), rp), this%Xh%zg(:,3), lz-1, 0, ht)
226 
227  ! And interpolate!
228  call triple_tensor_product(tmp(:,1), x, y, z, lx, hr, hs, ht)
229 
230  if (n .eq. 1) then
231  res(1)%x = tmp(:, 1)
232  return
233  end if
234 
235 
236  !
237  ! Loop through the rest of the points
238  !
239  do i = 2, n
240 
241  ! If the coordinates are different, then recompute weights
242  if ( .not. abscmp(rst(i)%x(1), rst(i-1)%x(1)) ) then
243  call fd_weights_full(real(rst(i)%x(1), rp), &
244  this%Xh%zg(:,1), lx-1, 0, hr)
245  end if
246  if ( .not. abscmp(rst(i)%x(2), rst(i-1)%x(2)) ) then
247  call fd_weights_full(real(rst(i)%x(2), rp), &
248  this%Xh%zg(:,2), ly-1, 0, hs)
249  end if
250  if ( .not. abscmp(rst(i)%x(3), rst(i-1)%x(3)) ) then
251  call fd_weights_full(real(rst(i)%x(3), rp), &
252  this%Xh%zg(:,3), lz-1, 0, ht)
253  end if
254 
255  ! And interpolate!
256  call triple_tensor_product(tmp(:,i), x, y, z, lx, hr, hs, ht)
257  end do
258 
259  ! Cast result to point_t dp
260  do i = 1, n
261  res(i)%x = tmp(:,i)
262  end do
263 
265 
274  function point_interpolator_interpolate_vector_jacobian(this, jac, rst, X, Y, Z) result(res)
275  class(point_interpolator_t), intent(in) :: this
276  real(kind=rp), intent(inout) :: jac(3,3)
277  type(point_t), intent(in) :: rst
278  real(kind=rp), intent(inout) :: x(this%Xh%lx, this%Xh%ly, this%Xh%lz)
279  real(kind=rp), intent(inout) :: y(this%Xh%lx, this%Xh%ly, this%Xh%lz)
280  real(kind=rp), intent(inout) :: z(this%Xh%lx, this%Xh%ly, this%Xh%lz)
281 
282  real(kind=rp) :: hr(this%Xh%lx, 2), hs(this%Xh%ly, 2), ht(this%Xh%lz, 2)
283  type(point_t) :: res
284  real(kind=rp) :: tmp(3)
285 
286  integer :: lx,ly,lz, i
287  lx = this%Xh%lx
288  ly = this%Xh%ly
289  lz = this%Xh%lz
290 
291  !
292  ! Compute weights
293  !
294  call fd_weights_full(real(rst%x(1), rp), this%Xh%zg(:,1), lx-1, 1, hr)
295  call fd_weights_full(real(rst%x(2), rp), this%Xh%zg(:,2), ly-1, 1, hs)
296  call fd_weights_full(real(rst%x(3), rp),this%Xh%zg(:,3), lz-1, 1, ht)
297 
298  !
299  ! Interpolate
300  !
301  call triple_tensor_product(tmp, x, y, z, lx, hr(:,1), hs(:,1), ht(:,1))
302  res%x = dble(tmp)! Cast from rp -> point_t dp
303 
304  !
305  ! Build jacobian
306  !
307 
308  ! d(x,y,z)/dr
309  call triple_tensor_product(tmp, x,y,z, lx, hr(:,2), hs(:,1), ht(:,1))
310  jac(1,:) = tmp
311 
312  ! d(x,y,z)/ds
313  call triple_tensor_product(tmp, x,y,z, lx, hr(:,1), hs(:,2), ht(:,1))
314  jac(2,:) = tmp
315 
316  ! d(x,y,z)/dt
317  call triple_tensor_product(tmp, x,y,z, lx, hr(:,1), hs(:,1), ht(:,2))
318  jac(3,:) = tmp
319 
321 
328  function point_interpolator_interpolate_jacobian(this, rst, X,Y,Z) result(jac)
329  class(point_interpolator_t), intent(in) :: this
330  type(point_t), intent(in) :: rst
331  real(kind=rp), intent(inout) :: x(this%Xh%lx, this%Xh%ly, this%Xh%lz)
332  real(kind=rp), intent(inout) :: y(this%Xh%lx, this%Xh%ly, this%Xh%lz)
333  real(kind=rp), intent(inout) :: z(this%Xh%lx, this%Xh%ly, this%Xh%lz)
334 
335  real(kind=rp) :: jac(3,3)
336  real(kind=rp) :: tmp(3)
337 
338  real(kind=rp) :: hr(this%Xh%lx, 2), hs(this%Xh%ly, 2), ht(this%Xh%lz, 2)
339  integer :: lx, ly, lz
340  lx = this%Xh%lx
341  ly = this%Xh%ly
342  lz = this%Xh%lz
343 
344  ! Weights
345  call fd_weights_full(real(rst%x(1), rp), this%Xh%zg(:,1), lx-1, 1, hr)
346  call fd_weights_full(real(rst%x(2), rp), this%Xh%zg(:,2), ly-1, 1, hs)
347  call fd_weights_full(real(rst%x(3), rp), this%Xh%zg(:,3), lz-1, 1, ht)
348 
349  ! d(x,y,z)/dr
350  call triple_tensor_product(tmp, x, y, z, lx, hr(:,2), hs(:,1), ht(:,1))
351  jac(1,:) = tmp
352 
353  ! d(x,y,z)/ds
354  call triple_tensor_product(tmp, x, y, z, lx, hr(:,1), hs(:,2), ht(:,1))
355  jac(2,:) = tmp
356 
357  ! d(x,y,z)/dt
358  call triple_tensor_product(tmp, x, y, z, lx, hr(:,1), hs(:,1), ht(:,2))
359  jac(3,:) = tmp
360 
362 
363 end module point_interpolator
double real
Definition: device_config.h:12
subroutine, public device_rzero(a_d, n)
Zero a real vector.
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
Fast diagonalization methods from NEKTON.
Definition: fast3d.f90:61
subroutine, public fd_weights_full(xi, x, n, m, c)
Compute finite-difference stencil weights for evaluating derivatives up to order at a point.
Definition: fast3d.f90:105
Definition: math.f90:60
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Routines to interpolate fields on a given element on a point in that element with given r,...
type(point_t) function point_interpolator_interpolate_vector_jacobian(this, jac, rst, X, Y, Z)
Interpolates a vector field and constructs the Jacobian at a point . Returns a vector .
subroutine point_interpolator_free(this)
Free pointers.
type(point_t) function, dimension(:), allocatable point_interpolator_interpolate_vector(this, rst, X, Y, Z)
Interpolates a vector field on a set of points . Returns an array of N points .
real(kind=rp) function, dimension(:), allocatable point_interpolator_interpolate_scalar(this, rst, X)
Interpolates a scalar field on a set of points . Returns a vector of N coordinates .
real(kind=rp) function, dimension(3, 3) point_interpolator_interpolate_jacobian(this, rst, X, Y, Z)
Constructs the Jacobian, returns a 3-by-3 array where .
subroutine point_interpolator_compute_weights(this, r, s, t, wr, ws, wt)
Computes interpolation weights for a list of points.
subroutine point_interpolator_init(this, Xh)
Initialization of point interpolation.
Implements a point.
Definition: point.f90:35
Defines a function space.
Definition: space.f90:34
integer, parameter, public gll
Definition: space.f90:48
integer, parameter, public gl
Definition: space.f90:48
Tensor operations.
Definition: tensor.f90:61
Utilities.
Definition: utils.f90:35
A point in with coordinates .
Definition: point.f90:43
Field interpolator to arbitrary points within an element. Tailored for experimentation,...
The function space for the SEM solution fields.
Definition: space.f90:62