Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
local_interpolation.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!
37 use space, only: space_t, gl, gll
38 use num_types, only: rp, xp
39 use point, only: point_t
40 use math, only: abscmp, neko_eps
41 use speclib
43 use utils, only: neko_error
44 use field, only: field_t
45 use field_list, only: field_list_t
46 use device
47 use math, only : matinv3, matinv39
48 use tensor_cpu
49 use device_math, only: device_rzero
51 use, intrinsic :: iso_c_binding
52 implicit none
53 private
54
60 type, public :: local_interpolator_t
62 type(space_t), pointer :: xh => null()
64 integer :: n_points
66 real(kind=rp), allocatable :: weights_r(:,:)
67 real(kind=rp), allocatable :: weights_s(:,:)
68 real(kind=rp), allocatable :: weights_t(:,:)
69 type(c_ptr) :: weights_r_d = c_null_ptr
70 type(c_ptr) :: weights_s_d = c_null_ptr
71 type(c_ptr) :: weights_t_d = c_null_ptr
72 contains
74 procedure, pass(this) :: init_3arrays => local_interpolator_init_3arrays
75 procedure, pass(this) :: init_1array => local_interpolator_init_1array
77 procedure, pass(this) :: free => local_interpolator_free
79 procedure, pass(this) :: evaluate => local_interpolator_evaluate
81 procedure, pass(this) :: compute_weights => local_interpolator_compute_weights
82 generic :: init => init_3arrays, init_1array
83
85
86contains
87
90 subroutine local_interpolator_init_3arrays(this, Xh, r, s, t, n_points)
91 class(local_interpolator_t), intent(inout), target :: this
92 type(space_t), intent(in), target :: Xh
93 integer, intent(in) :: n_points
94 real(kind=rp) :: r(n_points), s(n_points), t(n_points)
95 integer :: size_weights
96 call this%free()
97 if ((xh%t .eq. gl) .or. (xh%t .eq. gll)) then
98 else
99 call neko_error('Unsupported interpolation')
100 end if
101
102 this%Xh => xh
103 this%n_points = n_points
104 allocate(this%weights_r(xh%lx,n_points))
105 allocate(this%weights_s(xh%ly,n_points))
106 allocate(this%weights_t(xh%lz,n_points))
107 call this%compute_weights(r, s, t)
108 size_weights = xh%lx * n_points
109
110 if (neko_bcknd_device .eq. 1) then
111 call device_map(this%weights_r, this%weights_r_d, size_weights)
112 call device_map(this%weights_s, this%weights_s_d, size_weights)
113 call device_map(this%weights_t, this%weights_t_d, size_weights)
114 call device_memcpy(this%weights_r, this%weights_r_d,&
115 size_weights, host_to_device, sync = .true.)
116 call device_memcpy(this%weights_s, this%weights_s_d,&
117 size_weights, host_to_device, sync = .true.)
118 call device_memcpy(this%weights_t, this%weights_t_d,&
119 size_weights, host_to_device, sync = .true.)
120 end if
121
123
126 subroutine local_interpolator_init_1array(this, Xh, rst, n_points)
127 class(local_interpolator_t), intent(inout), target :: this
128 type(space_t), intent(in), target :: Xh
129 integer, intent(in) :: n_points
130 real(kind=rp), intent(in) :: rst(3,n_points)
131 real(kind=rp), allocatable :: r(:), s(:), t(:)
132 integer :: i
133
134 if (allocated(r)) deallocate(r)
135 allocate(r(n_points))
136 if (allocated(s)) deallocate(s)
137 allocate(s(n_points))
138 if (allocated(t)) deallocate(t)
139 allocate(t(n_points))
140
141 do i = 1, n_points
142 r(i) = rst(1,i)
143 s(i) = rst(2,i)
144 t(i) = rst(3,i)
145 end do
146
147 call this%init_3arrays(xh, r, s, t, n_points)
148
149 deallocate(r,s,t)
150
151 end subroutine local_interpolator_init_1array
152
153
155 subroutine local_interpolator_free(this)
156 class(local_interpolator_t), intent(inout) :: this
157
158 if (associated(this%Xh)) this%Xh => null()
159
160 if(allocated(this%weights_r)) deallocate(this%weights_r)
161 if(allocated(this%weights_s)) deallocate(this%weights_s)
162 if(allocated(this%weights_t)) deallocate(this%weights_t)
163 if (c_associated(this%weights_r_d)) then
164 call device_free(this%weights_r_d)
165 end if
166 if (c_associated(this%weights_s_d)) then
167 call device_free(this%weights_s_d)
168 end if
169 if (c_associated(this%weights_t_d)) then
170 call device_free(this%weights_t_d)
171 end if
172 end subroutine local_interpolator_free
173
184 subroutine local_interpolator_compute_weights(this, r, s, t)
185 class(local_interpolator_t), intent(inout) :: this
186 real(kind=rp), intent(in) :: r(:), s(:), t(:)
187
188 integer :: N, i, lx
189 lx = this%Xh%lx
190 n = size(r)
191
192 do i = 1, n
193 if ((r(i) <= 1.1_rp .and. r(i) >= -1.1_rp) .and. &
194 (s(i) <= 1.1_rp .and. s(i) >= -1.1_rp) .and. &
195 (t(i) <= 1.1_rp .and. t(i) >= -1.1_rp)) then
196 call fd_weights_full(r(i), this%Xh%zg(:,1), lx-1, 0, this%weights_r(:,i))
197 call fd_weights_full(s(i), this%Xh%zg(:,2), lx-1, 0, this%weights_s(:,i))
198 call fd_weights_full(t(i), this%Xh%zg(:,3), lx-1, 0, this%weights_t(:,i))
199 else
200 this%weights_r(:,i) = 0.0_rp
201 this%weights_s(:,i) = 0.0_rp
202 this%weights_t(:,i) = 0.0_rp
203 end if
204
205 end do
206
208
221 subroutine local_interpolator_evaluate(this, interp_values, el_list, field, &
222 nel, on_host)
223 class(local_interpolator_t), intent(inout) :: this
224 integer, intent(in) :: el_list(this%n_points)
225 integer, intent(in) :: nel
226 real(kind=rp), intent(inout) :: interp_values(this%n_points)
227 real(kind=rp), intent(inout) :: field(this%Xh%lxyz, nel)
228 logical, intent(in) :: on_host
229
230 call tnsr3d_el_list(interp_values, 1, field, this%Xh%lx, &
231 this%weights_r, this%weights_s, this%weights_t, el_list, &
232 this%n_points, on_host)
233
234 end subroutine local_interpolator_evaluate
235
242 subroutine jacobian(jac, rst, x, y, z, n_pts, Xh)
243 integer, intent(in) :: n_pts
244 real(kind=rp), intent(inout) :: rst(3, n_pts)
245 type(space_t), intent(inout) :: xh
246 real(kind=rp), intent(inout) :: x(xh%lx, xh%ly, xh%lz, n_pts)
247 real(kind=rp), intent(inout) :: y(xh%lx, xh%ly, xh%lz, n_pts)
248 real(kind=rp), intent(inout) :: z(xh%lx, xh%ly, xh%lz, n_pts)
249 real(kind=rp), intent(out) :: jac(3,3, n_pts)
250 real(kind=rp) :: tmp(3)
251 real(kind=rp) :: hr(xh%lx, 2), hs(xh%ly, 2), ht(xh%lz, 2)
252 integer :: lx, ly, lz, i
253 lx = xh%lx
254 ly = xh%ly
255 lz = xh%lz
256
257 do i = 1, n_pts
258 ! Weights
259 call fd_weights_full(rst(1,i), xh%zg(:,1), lx-1, 1, hr)
260 call fd_weights_full(rst(2,i), xh%zg(:,2), ly-1, 1, hs)
261 call fd_weights_full(rst(3,i), xh%zg(:,3), lz-1, 1, ht)
262
263 ! d(x,y,z)/dr
264 call triple_tensor_product(tmp(1), x(:,:,:,i), lx, hr(:,2), hs(:,1), &
265 ht(:,1))
266 jac(1,1,i) = tmp(1)
267 call triple_tensor_product(tmp(1), y(:,:,:,i), lx, hr(:,2), hs(:,1), &
268 ht(:,1))
269 jac(1,2,i) = tmp(1)
270 call triple_tensor_product(tmp(1), z(:,:,:,i), lx, hr(:,2), hs(:,1), &
271 ht(:,1))
272 jac(1,3,i) = tmp(1)
273
274 ! d(x,y,z)/ds
275 call triple_tensor_product(tmp, x(:,:,:,i), y(:,:,:,i), z(:,:,:,i), lx, &
276 hr(:,1), hs(:,2), ht(:,1))
277 jac(2,:,i) = tmp
278
279 ! d(x,y,z)/dt
280 call triple_tensor_product(tmp, x(:,:,:,i), y(:,:,:,i), z(:,:,:,i), lx, &
281 hr(:,1), hs(:,1), ht(:,2))
282 jac(3,:,i) = tmp
283 end do
284
285 end subroutine jacobian
286
287 subroutine jacobian_inverse(jacinv, rst, x, y, z, n_pts, Xh)
288 integer :: n_pts
289 real(kind=rp), intent(inout) :: rst(3, n_pts)
290 type(space_t), intent(inout) :: xh
291 real(kind=rp), intent(inout) :: x(xh%lx, xh%ly, xh%lz, n_pts)
292 real(kind=rp), intent(inout) :: y(xh%lx, xh%ly, xh%lz, n_pts)
293 real(kind=rp), intent(inout) :: z(xh%lx, xh%ly, xh%lz, n_pts)
294 real(kind=rp), intent(out) :: jacinv(3,3, n_pts)
295 real(kind=rp) :: tmp(3,3)
296 integer :: i
297
298 call jacobian(jacinv, rst, x, y, z, n_pts, xh)
299
300 do i = 1, n_pts
301 tmp = matinv3(real(jacinv(:,:,3),xp))
302 jacinv(:,:,i) = tmp
303 end do
304
305 end subroutine jacobian_inverse
306
307end module local_interpolation
double real
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
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
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
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
subroutine, public setup_intp(jh, jht, z_to, z_from, n_to, n_from, derivative)
Compute interpolation weights for points z_to using values at points z_from.
Definition fast3d.f90:244
Defines a field.
Definition field.f90:34
Routines to obtain interpolated values on a set of points with known rst coordinates in elements loca...
subroutine local_interpolator_compute_weights(this, r, s, t)
Computes interpolation weights for a list of points.
subroutine local_interpolator_free(this)
Free pointers.
subroutine jacobian(jac, rst, x, y, z, n_pts, xh)
Constructs the Jacobian, returns a 3-by-3 times number of points where .
subroutine local_interpolator_init_3arrays(this, xh, r, s, t, n_points)
Initialization of point interpolation.
subroutine local_interpolator_init_1array(this, xh, rst, n_points)
Initialization of point interpolation.
subroutine jacobian_inverse(jacinv, rst, x, y, z, n_pts, xh)
subroutine local_interpolator_evaluate(this, interp_values, el_list, field, nel, on_host)
Interpolates a list of fields based on a set of element ids.
Definition math.f90:60
real(rp) function, dimension(3, 3), public matinv39(a11, a12, a13, a21, a22, a23, a31, a32, a33)
Definition math.f90:1486
real(kind=xp) function, dimension(3, 3), public matinv3(a)
Performs a direct calculation of the inverse of a 3×3 matrix. M33INV and M44INV by David G....
Definition math.f90:1506
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:69
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public xp
Definition num_types.f90:14
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements a point.
Definition point.f90:35
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:49
integer, parameter, public gl
Definition space.f90:49
LIBRARY ROUTINES FOR SPECTRAL METHODS.
Definition speclib.f90:148
Tensor operations.
Definition tensor.f90:61
subroutine, public tnsr3d_el_list(v, nv, u, nu, a, bt, ct, el_list, n_pt, on_host)
Tensor product performed on a subset of the elements.
Definition tensor.f90:191
subroutine, public tnsr3d(v, nv, u, nu, a, bt, ct, nelv)
Tensor product performed on nelv elements.
Definition tensor.f90:234
Utilities.
Definition utils.f90:35
field_list_t, To be able to group fields together
Interpolation on a set of points with known rst coordinates in elements local to this process....
A point in with coordinates .
Definition point.f90:43
The function space for the SEM solution fields.
Definition space.f90:63