88 type(
space_t),
intent(in),
target :: Xh
90 if ((xh%t .eq.
gl) .or. (xh%t .eq.
gll))
then
103 if (
associated(this%Xh)) this%Xh => null()
119 real(kind=
rp),
intent(in) :: r(:), s(:), t(:)
120 real(kind=
rp),
intent(inout) :: wr(:,:), ws(:,:), wt(:,:)
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(:)
145 real(kind=
rp) :: hr(this%Xh%lx), hs(this%Xh%ly), ht(this%Xh%lz)
146 integer :: lx,ly,lz, i
173 if ( .not.
abscmp(rst(i)%x(1), rst(i-1)%x(1)) )
then
175 this%Xh%zg(:,1), lx-1, 0, hr)
177 if ( .not.
abscmp(rst(i)%x(2), rst(i-1)%x(2)) )
then
179 this%Xh%zg(:,2), ly-1, 0, hs)
181 if ( .not.
abscmp(rst(i)%x(3), rst(i-1)%x(3)) )
then
183 this%Xh%zg(:,3), lz-1, 0, ht)
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)
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
242 if ( .not.
abscmp(rst(i)%x(1), rst(i-1)%x(1)) )
then
244 this%Xh%zg(:,1), lx-1, 0, hr)
246 if ( .not.
abscmp(rst(i)%x(2), rst(i-1)%x(2)) )
then
248 this%Xh%zg(:,2), ly-1, 0, hs)
250 if ( .not.
abscmp(rst(i)%x(3), rst(i-1)%x(3)) )
then
252 this%Xh%zg(:,3), lz-1, 0, ht)
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)
282 real(kind=
rp) :: hr(this%Xh%lx, 2), hs(this%Xh%ly, 2), ht(this%Xh%lz, 2)
284 real(kind=
rp) :: tmp(3)
286 integer :: lx,ly,lz, i
301 call triple_tensor_product(tmp, x, y, z, lx, hr(:,1), hs(:,1), ht(:,1))
309 call triple_tensor_product(tmp, x,y,z, lx, hr(:,2), hs(:,1), ht(:,1))
313 call triple_tensor_product(tmp, x,y,z, lx, hr(:,1), hs(:,2), ht(:,1))
317 call triple_tensor_product(tmp, x,y,z, lx, hr(:,1), hs(:,1), ht(:,2))
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)
335 real(kind=
rp) :: jac(3,3)
336 real(kind=
rp) :: tmp(3)
338 real(kind=
rp) :: hr(this%Xh%lx, 2), hs(this%Xh%ly, 2), ht(this%Xh%lz, 2)
339 integer :: lx, ly, lz
350 call triple_tensor_product(tmp, x, y, z, lx, hr(:,2), hs(:,1), ht(:,1))
354 call triple_tensor_product(tmp, x, y, z, lx, hr(:,1), hs(:,2), ht(:,1))
358 call triple_tensor_product(tmp, x, y, z, lx, hr(:,1), hs(:,1), ht(:,2))
subroutine, public device_rzero(a_d, n)
Zero a real vector.
Device abstraction, common interface for various accelerators.
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.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Routines to interpolate fields on a given element on a point in that element with given r,...
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 .
subroutine point_interpolator_free(this)
Free pointers.
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 .
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_init(this, xh)
Initialization of point interpolation.
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 .
subroutine point_interpolator_compute_weights(this, r, s, t, wr, ws, wt)
Computes interpolation weights for a list of points.
Defines a function space.
integer, parameter, public gll
integer, parameter, public gl
A point in with coordinates .
Field interpolator to arbitrary points within an element. Tailored for experimentation,...
The function space for the SEM solution fields.