Neko 1.99.3
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
42 use utils, only : neko_error
43 use field, only : field_t
44 use field_list, only : field_list_t
45 use device
46 use math, only : matinv3, matinv39
47 use device_math, only : device_rzero
49 use, intrinsic :: iso_c_binding
50 implicit none
51 private
52
58 type, public :: local_interpolator_t
60 type(space_t), pointer :: xh => null()
62 integer :: n_points
64 real(kind=rp), allocatable :: weights_r(:,:)
65 real(kind=rp), allocatable :: weights_s(:,:)
66 real(kind=rp), allocatable :: weights_t(:,:)
67 type(c_ptr) :: weights_r_d = c_null_ptr
68 type(c_ptr) :: weights_s_d = c_null_ptr
69 type(c_ptr) :: weights_t_d = c_null_ptr
70 contains
72 procedure, pass(this) :: init_3arrays => local_interpolator_init_3arrays
73 procedure, pass(this) :: init_1array => local_interpolator_init_1array
75 procedure, pass(this) :: free => local_interpolator_free
77 procedure, pass(this) :: evaluate => local_interpolator_evaluate
79 procedure, pass(this) :: compute_weights => &
81 generic :: init => init_3arrays, init_1array
82
84
85contains
86
89 subroutine local_interpolator_init_3arrays(this, Xh, r, s, t, n_points)
90 class(local_interpolator_t), intent(inout), target :: this
91 type(space_t), intent(in), target :: Xh
92 integer, intent(in) :: n_points
93 real(kind=rp) :: r(n_points), s(n_points), t(n_points)
94 integer :: size_weights
95 call this%free()
96 if ((xh%t .eq. gl) .or. (xh%t .eq. gll)) then
97 else
98 call neko_error('Unsupported interpolation')
99 end if
100
101 this%Xh => xh
102 this%n_points = n_points
103 allocate(this%weights_r(xh%lx, n_points))
104 allocate(this%weights_s(xh%ly, n_points))
105 allocate(this%weights_t(xh%lz, n_points))
106 call this%compute_weights(r, s, t)
107 size_weights = xh%lx * n_points
108
109 if (neko_bcknd_device .eq. 1) then
110 call device_map(this%weights_r, this%weights_r_d, size_weights)
111 call device_map(this%weights_s, this%weights_s_d, size_weights)
112 call device_map(this%weights_t, this%weights_t_d, size_weights)
113 call device_memcpy(this%weights_r, this%weights_r_d,&
114 size_weights, host_to_device, sync = .true.)
115 call device_memcpy(this%weights_s, this%weights_s_d,&
116 size_weights, host_to_device, sync = .true.)
117 call device_memcpy(this%weights_t, this%weights_t_d,&
118 size_weights, host_to_device, sync = .true.)
119 end if
120
122
125 subroutine local_interpolator_init_1array(this, Xh, rst, n_points)
126 class(local_interpolator_t), intent(inout), target :: this
127 type(space_t), intent(in), target :: Xh
128 integer, intent(in) :: n_points
129 real(kind=rp), intent(in) :: rst(3, n_points)
130 real(kind=rp), allocatable :: r(:), s(:), t(:)
131 integer :: i
132
133 if (allocated(r)) deallocate(r)
134 allocate(r(n_points))
135 if (allocated(s)) deallocate(s)
136 allocate(s(n_points))
137 if (allocated(t)) deallocate(t)
138 allocate(t(n_points))
139
140 do i = 1, n_points
141 r(i) = rst(1,i)
142 s(i) = rst(2,i)
143 t(i) = rst(3,i)
144 end do
145
146 call this%init_3arrays(xh, r, s, t, n_points)
147
148 deallocate(r,s,t)
149
150 end subroutine local_interpolator_init_1array
151
152
154 subroutine local_interpolator_free(this)
155 class(local_interpolator_t), intent(inout) :: this
156
157 if (associated(this%Xh)) this%Xh => null()
158
159 if (allocated(this%weights_r)) then
160 if (neko_bcknd_device .eq. 1) then
161 call device_unmap(this%weights_r, this%weights_r_d)
162 end if
163 deallocate(this%weights_r)
164 end if
165 if (allocated(this%weights_s)) then
166 if (neko_bcknd_device .eq. 1) then
167 call device_unmap(this%weights_s, this%weights_s_d)
168 end if
169 deallocate(this%weights_s)
170 end if
171 if (allocated(this%weights_t)) then
172 if (neko_bcknd_device .eq. 1) then
173 call device_unmap(this%weights_t, this%weights_t_d)
174 end if
175 deallocate(this%weights_t)
176 end if
177 end subroutine local_interpolator_free
178
189 subroutine local_interpolator_compute_weights(this, r, s, t)
190 class(local_interpolator_t), intent(inout) :: this
191 real(kind=rp), intent(in) :: r(:), s(:), t(:)
192
193 integer :: N, i, lx
194 lx = this%Xh%lx
195 n = size(r)
196
197 do i = 1, n
198 if ((r(i) <= 1.1_rp .and. r(i) >= -1.1_rp) .and. &
199 (s(i) <= 1.1_rp .and. s(i) >= -1.1_rp) .and. &
200 (t(i) <= 1.1_rp .and. t(i) >= -1.1_rp)) then
201 call fd_weights_full(r(i), this%Xh%zg(:, 1), lx-1, 0, &
202 this%weights_r(:, i))
203 call fd_weights_full(s(i), this%Xh%zg(:, 2), lx-1, 0, &
204 this%weights_s(:, i))
205 call fd_weights_full(t(i), this%Xh%zg(:, 3), lx-1, 0, &
206 this%weights_t(:, i))
207 else
208 this%weights_r(:,i) = 0.0_rp
209 this%weights_s(:,i) = 0.0_rp
210 this%weights_t(:,i) = 0.0_rp
211 end if
212
213 end do
214
216
229 subroutine local_interpolator_evaluate(this, interp_values, el_list, field, &
230 nel, on_host)
231 class(local_interpolator_t), intent(inout) :: this
232 integer, intent(in) :: el_list(this%n_points)
233 integer, intent(in) :: nel
234 real(kind=rp), intent(inout) :: interp_values(this%n_points)
235 real(kind=rp), intent(inout) :: field(this%Xh%lxyz, nel)
236 logical, intent(in) :: on_host
237
238 call tnsr3d_el_list(interp_values, 1, field, this%Xh%lx, &
239 this%weights_r, this%weights_s, this%weights_t, el_list, &
240 this%n_points, on_host)
241
242 end subroutine local_interpolator_evaluate
243
250 subroutine jacobian(jac, rst, x, y, z, n_pts, Xh)
251 integer, intent(in) :: n_pts
252 real(kind=rp), intent(inout) :: rst(3, n_pts)
253 type(space_t), intent(inout) :: xh
254 real(kind=rp), intent(inout) :: x(xh%lx, xh%ly, xh%lz, n_pts)
255 real(kind=rp), intent(inout) :: y(xh%lx, xh%ly, xh%lz, n_pts)
256 real(kind=rp), intent(inout) :: z(xh%lx, xh%ly, xh%lz, n_pts)
257 real(kind=rp), intent(out) :: jac(3,3, n_pts)
258 real(kind=rp) :: tmp(3)
259 real(kind=rp) :: hr(xh%lx, 2), hs(xh%ly, 2), ht(xh%lz, 2)
260 integer :: lx, ly, lz, i
261 lx = xh%lx
262 ly = xh%ly
263 lz = xh%lz
264
265 do i = 1, n_pts
266 ! Weights
267 call fd_weights_full(rst(1,i), xh%zg(:,1), lx-1, 1, hr)
268 call fd_weights_full(rst(2,i), xh%zg(:,2), ly-1, 1, hs)
269 call fd_weights_full(rst(3,i), xh%zg(:,3), lz-1, 1, ht)
270
271 ! d(x,y,z)/dr
272 call triple_tensor_product(tmp(1), x(:,:,:,i), lx, hr(:,2), hs(:,1), &
273 ht(:,1))
274 jac(1,1,i) = tmp(1)
275 call triple_tensor_product(tmp(1), y(:,:,:,i), lx, hr(:,2), hs(:,1), &
276 ht(:,1))
277 jac(1,2,i) = tmp(1)
278 call triple_tensor_product(tmp(1), z(:,:,:,i), lx, hr(:,2), hs(:,1), &
279 ht(:,1))
280 jac(1,3,i) = tmp(1)
281
282 ! d(x,y,z)/ds
283 call triple_tensor_product(tmp, x(:,:,:,i), y(:,:,:,i), z(:,:,:,i), lx, &
284 hr(:,1), hs(:,2), ht(:,1))
285 jac(2,:,i) = tmp
286
287 ! d(x,y,z)/dt
288 call triple_tensor_product(tmp, x(:,:,:,i), y(:,:,:,i), z(:,:,:,i), lx, &
289 hr(:,1), hs(:,1), ht(:,2))
290 jac(3,:,i) = tmp
291 end do
292
293 end subroutine jacobian
294
295 subroutine jacobian_inverse(jacinv, rst, x, y, z, n_pts, Xh)
296 integer :: n_pts
297 real(kind=rp), intent(inout) :: rst(3, n_pts)
298 type(space_t), intent(inout) :: xh
299 real(kind=rp), intent(inout) :: x(xh%lx, xh%ly, xh%lz, n_pts)
300 real(kind=rp), intent(inout) :: y(xh%lx, xh%ly, xh%lz, n_pts)
301 real(kind=rp), intent(inout) :: z(xh%lx, xh%ly, xh%lz, n_pts)
302 real(kind=rp), intent(out) :: jacinv(3,3, n_pts)
303 real(kind=rp) :: tmp(3,3)
304 integer :: i
305
306 call jacobian(jacinv, rst, x, y, z, n_pts, xh)
307
308 do i = 1, n_pts
309 tmp = matinv3(real(jacinv(:, :, 3), xp))
310 jacinv(:, :, i) = tmp
311 end do
312
313 end subroutine jacobian_inverse
314
315end module local_interpolation
double real
Map a Fortran array to a device (allocate and associate)
Definition device.F90:78
Copy data between host and device (or device and device)
Definition device.F90:72
Unmap a Fortran array from a device (deassociate and free)
Definition device.F90:84
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:48
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:1780
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:1800
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:70
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
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