Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
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
39 use point, only: point_t
40 use math, only: abscmp
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 device_math, only: device_rzero
48 use, intrinsic :: iso_c_binding
49 implicit none
50 private
51
57 type, public :: local_interpolator_t
59 type(space_t), pointer :: xh => null()
61 integer :: n_points
63 real(kind=rp), allocatable :: weights_r(:,:)
64 real(kind=rp), allocatable :: weights_s(:,:)
65 real(kind=rp), allocatable :: weights_t(:,:)
66 type(c_ptr) :: weights_r_d = c_null_ptr
67 type(c_ptr) :: weights_s_d = c_null_ptr
68 type(c_ptr) :: weights_t_d = c_null_ptr
69 contains
71 procedure, pass(this) :: init => local_interpolator_init
73 procedure, pass(this) :: free => local_interpolator_free
75 procedure, pass(this) :: evaluate => local_interpolator_evaluate
77 procedure, pass(this) :: compute_weights => local_interpolator_compute_weights
78
80
81contains
82
85 subroutine local_interpolator_init(this, Xh, r, s, t, n_points)
86 class(local_interpolator_t), intent(inout), target :: this
87 type(space_t), intent(in), target :: Xh
88 integer, intent(in) :: n_points
89 real(kind=rp) :: r(n_points), s(n_points), t(n_points)
90 integer :: size_weights
91 call this%free()
92 if ((xh%t .eq. gl) .or. (xh%t .eq. gll)) then
93 else
94 call neko_error('Unsupported interpolation')
95 end if
96
97 this%Xh => xh
98 this%n_points = n_points
99 allocate(this%weights_r(xh%lx,n_points))
100 allocate(this%weights_s(xh%ly,n_points))
101 allocate(this%weights_t(xh%lz,n_points))
102 call this%compute_weights(r, s, t)
103 size_weights = xh%lx * n_points
104
105 if (neko_bcknd_device .eq. 1) then
106 call device_map(this%weights_r, this%weights_r_d, size_weights)
107 call device_map(this%weights_s, this%weights_s_d, size_weights)
108 call device_map(this%weights_t, this%weights_t_d, size_weights)
109 call device_memcpy(this%weights_r, this%weights_r_d,&
110 size_weights, host_to_device, sync = .true.)
111 call device_memcpy(this%weights_s, this%weights_s_d,&
112 size_weights, host_to_device, sync = .true.)
113 call device_memcpy(this%weights_t, this%weights_t_d,&
114 size_weights, host_to_device, sync = .true.)
115 end if
116
117 end subroutine local_interpolator_init
118
120 subroutine local_interpolator_free(this)
121 class(local_interpolator_t), intent(inout) :: this
122
123 if (associated(this%Xh)) this%Xh => null()
124
125 if(allocated(this%weights_r)) deallocate(this%weights_r)
126 if(allocated(this%weights_s)) deallocate(this%weights_s)
127 if(allocated(this%weights_t)) deallocate(this%weights_t)
128 if (c_associated(this%weights_r_d)) then
129 call device_free(this%weights_r_d)
130 end if
131 if (c_associated(this%weights_s_d)) then
132 call device_free(this%weights_s_d)
133 end if
134 if (c_associated(this%weights_t_d)) then
135 call device_free(this%weights_t_d)
136 end if
137 end subroutine local_interpolator_free
138
149 subroutine local_interpolator_compute_weights(this, r, s, t)
150 class(local_interpolator_t), intent(inout) :: this
151 real(kind=rp), intent(in) :: r(:), s(:), t(:)
152
153 integer :: N, i, lx
154 lx = this%Xh%lx
155 n = size(r)
156
157 do i = 1, n
158 call fd_weights_full(r(i), this%Xh%zg(:,1), lx-1, 0, this%weights_r(:,i))
159 call fd_weights_full(s(i), this%Xh%zg(:,2), lx-1, 0, this%weights_s(:,i))
160 call fd_weights_full(t(i), this%Xh%zg(:,3), lx-1, 0, this%weights_t(:,i))
161 end do
162
164
177 subroutine local_interpolator_evaluate(this, interp_values, el_list, field,nel)
178 class(local_interpolator_t), intent(inout) :: this
179 integer, intent(in) :: el_list(this%n_points)
180 integer, intent(in) :: nel
181 real(kind=rp), intent(inout) :: interp_values(this%n_points)
182 real(kind=rp), intent(inout) :: field(this%Xh%lxyz, nel)
183
184
185 call tnsr3d_el_list(interp_values, 1, field, this%Xh%lx, &
186 this%weights_r, this%weights_s, this%weights_t, el_list, this%n_points)
187
188 end subroutine local_interpolator_evaluate
189
190
191end module local_interpolation
Map a Fortran array to a device (allocate and associate)
Definition device.F90:71
Copy data between host and device (or device and device)
Definition device.F90:65
subroutine, public device_rzero(a_d, n)
Zero a real vector.
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:46
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:200
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:243
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_init(this, xh, r, s, t, n_points)
Initialization of point interpolation.
subroutine local_interpolator_evaluate(this, interp_values, el_list, field, nel)
Interpolates a list of fields based on a set of element ids.
subroutine local_interpolator_free(this)
Free pointers.
Definition math.f90:60
Build configurations.
integer, parameter neko_bcknd_device
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:48
integer, parameter, public gl
Definition space.f90:48
Tensor operations.
Definition tensor.f90:61
subroutine, public tnsr3d_el_list(v, nv, u, nu, a, bt, ct, el_list, n_pt)
Tensor product performed on a subset of the elements.
Definition tensor.f90:188
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:62