49 use,
intrinsic :: iso_c_binding
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
79 procedure, pass(this) :: compute_weights => &
81 generic :: init => init_3arrays, init_1array
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
96 if ((xh%t .eq.
gl) .or. (xh%t .eq.
gll))
then
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
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)
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(:)
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))
146 call this%init_3arrays(xh, r, s, t, n_points)
157 if (
associated(this%Xh)) this%Xh => null()
159 if (
allocated(this%weights_r))
then
163 deallocate(this%weights_r)
165 if (
allocated(this%weights_s))
then
169 deallocate(this%weights_s)
171 if (
allocated(this%weights_t))
then
175 deallocate(this%weights_t)
191 real(kind=
rp),
intent(in) :: r(:), s(:), t(:)
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
202 this%weights_r(:, i))
204 this%weights_s(:, i))
206 this%weights_t(:, i))
208 this%weights_r(:,i) = 0.0_rp
209 this%weights_s(:,i) = 0.0_rp
210 this%weights_t(:,i) = 0.0_rp
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
239 this%weights_r, this%weights_s, this%weights_t, el_list, &
240 this%n_points, on_host)
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
272 call triple_tensor_product(tmp(1), x(:,:,:,i), lx, hr(:,2), hs(:,1), &
275 call triple_tensor_product(tmp(1), y(:,:,:,i), lx, hr(:,2), hs(:,1), &
278 call triple_tensor_product(tmp(1), z(:,:,:,i), lx, hr(:,2), hs(:,1), &
283 call triple_tensor_product(tmp, x(:,:,:,i), y(:,:,:,i), z(:,:,:,i), lx, &
284 hr(:,1), hs(:,2), ht(:,1))
288 call triple_tensor_product(tmp, x(:,:,:,i), y(:,:,:,i), z(:,:,:,i), lx, &
289 hr(:,1), hs(:,1), ht(:,2))
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)
306 call jacobian(jacinv, rst, x, y, z, n_pts, xh)
310 jacinv(:, :, i) = tmp
Map a Fortran array to a device (allocate and associate)
Copy data between host and device (or device and device)
Unmap a Fortran array from a device (deassociate and free)
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
Fast diagonalization methods from NEKTON.
subroutine, public fd_weights_full(xi, x, n, m, c)
Compute finite-difference stencil weights for evaluating derivatives up to order at a point.
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.
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.
real(rp) function, dimension(3, 3), public matinv39(a11, a12, a13, a21, a22, a23, a31, a32, a33)
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....
real(kind=rp), parameter, public neko_eps
Machine epsilon .
integer, parameter neko_bcknd_device
integer, parameter, public xp
integer, parameter, public rp
Global precision used in computations.
Defines a function space.
integer, parameter, public gll
integer, parameter, public gl
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.
subroutine, public tnsr3d(v, nv, u, nu, a, bt, ct, nelv)
Tensor product performed on nelv elements.
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 .
The function space for the SEM solution fields.