73 use,
intrinsic :: iso_c_binding, only : c_ptr
93 subroutine tensr3(v, nv, u, nu, A, Bt, Ct, w)
96 real(kind=
rp),
intent(inout) :: v(nv, nv, nv)
97 real(kind=
rp),
intent(inout) :: u(nu, nu, nu)
98 real(kind=
rp),
intent(inout) :: w(nu*nu*nv)
99 real(kind=
rp),
intent(inout) :: a(nv, nu)
100 real(kind=
rp),
intent(inout) :: bt(nu, nv)
101 real(kind=
rp),
intent(inout) :: ct(nu, nv)
102 integer :: j, k, l, nunu, nvnv, nunv
110 call mxm(a, nv, u, nu, v, nunu)
114 call mxm(v(k, 1, 1), nv, bt, nu, w(l), nv)
118 call mxm(w, nvnv, ct, nu, v, nv)
123 subroutine trsp(a, lda, b, ldb)
124 integer,
intent(in) :: lda
125 integer,
intent(in) :: ldb
126 real(kind=
rp),
intent(inout) :: a(lda, ldb)
127 real(kind=
rp),
intent(in) :: b(ldb, lda)
140 integer,
intent(in) :: n
141 real(kind=
rp),
intent(inout) :: a(n, n)
157 integer,
intent(in) :: nv, nu
158 real(kind=
rp),
intent(inout) :: v(nv*nv), u(nu*nu)
159 real(kind=
rp),
intent(inout) :: a(nv, nu), bt(nu, nv)
174 integer,
intent(in) :: nv, nu
175 real(kind=
rp),
intent(inout) :: v(nv*nv*nv), u(nu*nu*nu)
176 real(kind=
rp),
intent(inout) :: a(nv,nu),bt(nu, nv),ct(nu,nv)
191 integer,
intent(in) :: nv, nu, n_pt, el_list(n_pt)
192 real(kind=
rp),
intent(inout) :: v(nv*nv*nv, n_pt), u(nu*nu*nu, 1)
193 real(kind=
rp),
intent(inout) :: a(nv, nu, n_pt)
194 real(kind=
rp),
intent(inout) :: bt(nu, nv, n_pt)
195 real(kind=
rp),
intent(inout) :: ct(nu, nv, n_pt)
196 logical,
intent(in) :: on_host
197 type(c_ptr) :: v_d, u_d, a_d, bt_d, ct_d, el_list_d
200 if (n_pt .eq. 0)
return
205 nu, a(1,1,i), bt(1,1,i), ct(1,1,i))
210 nu, a(1,1,i), bt(1,1,i), ct(1,1,i))
220 a_d, bt_d, ct_d, el_list_d, n_pt)
225 nu, a(1,1,i), bt(1,1,i), ct(1,1,i))
233 subroutine tnsr3d(v, nv, u, nu, A, Bt, Ct, nelv)
234 integer,
intent(in) :: nv, nu, nelv
235 real(kind=
rp),
intent(inout) :: v(nv*nv*nv,nelv)
236 real(kind=
rp),
intent(in) :: u(nu*nu*nu,nelv)
237 real(kind=
rp),
intent(in) :: a(nv,nu), bt(nu, nv), ct(nu,nv)
238 type(c_ptr) :: v_d, u_d, a_d, bt_d, ct_d
240 if (nelv .eq. 0)
return
243 call tnsr3d_sx(v, nv, u, nu, a, bt, ct, nelv)
254 call tnsr3d_cpu(v, nv, u, nu, a, bt, ct, nelv)
261 integer,
intent(inout) :: nv, nu, nelv
262 real(kind=
rp),
intent(inout) :: v(nv*nv*nv*nelv)
263 real(kind=
rp),
intent(inout) :: a(nv,nu), bt(nu, nv), ct(nu,nv)
265 if (nelv .eq. 0)
return
281 integer,
intent(in) :: nx, ny, nz
282 real(kind=
rp),
intent(in) :: h1(nx), h2(ny), h3(nz)
283 real(kind=
rp),
intent(inout) :: s(nx, ny, nz)
285 integer :: ix, iy, iz
291 s(ix,iy,iz) = s(ix,iy,iz)+hh*h1(ix)
313 real(kind=
rp),
intent(inout) :: v
314 integer,
intent(in) :: nu
315 real(kind=
rp),
intent(inout) :: u(nu, nu, nu)
316 real(kind=
rp),
intent(inout) :: hr(nu)
317 real(kind=
rp),
intent(inout) :: hs(nu)
318 real(kind=
rp),
intent(inout) :: ht(nu)
322 real(kind=
rp) :: vv(1)
349 real(kind=
rp),
intent(inout) :: v(3)
350 integer,
intent(in) :: nu
351 real(kind=
rp),
intent(inout) :: u1(nu, nu, nu)
352 real(kind=
rp),
intent(inout) :: u2(nu, nu, nu)
353 real(kind=
rp),
intent(inout) :: u3(nu, nu, nu)
354 real(kind=
rp),
intent(inout) :: hr(nu)
355 real(kind=
rp),
intent(inout) :: hs(nu)
356 real(kind=
rp),
intent(inout) :: ht(nu)
Return the device pointer for an associated Fortran array.
Device abstraction, common interface for various accelerators.
Wrapper for all matrix-matrix product implementations.
subroutine, public mxm(a, n1, b, n2, c, n3)
Compute matrix-matrix product for contiguously packed matrices A,B, and C.
integer, parameter neko_bcknd_sx
integer, parameter neko_bcknd_device
integer, parameter neko_bcknd_xsmm
integer, parameter, public rp
Global precision used in computations.
subroutine, public tnsr3d_cpu(v, nv, u, nu, a, bt, ct, nelv)
subroutine, public tnsr2d_el_cpu(v, nv, u, nu, a, bt)
subroutine, public tnsr1_3d_cpu(v, nv, nu, a, bt, ct, nelv)
subroutine, public tnsr3d_el_cpu(v, nv, u, nu, a, bt, ct)
subroutine, public tnsr3d_el_list_device(v_d, nv, u_d, nu, a_d, bt_d, ct_d, elements, n_points)
subroutine, public tnsr3d_device(v_d, nv, u_d, nu, a_d, bt_d, ct_d, nelv)
Tensor operations SX-Aurora backend.
subroutine, public tnsr1_3d_sx(v, nv, nu, a, bt, ct, nelv)
subroutine, public tnsr3d_sx(v, nv, u, nu, a, bt, ct, nelv)
subroutine, public tnsr2d_el_sx(v, nv, u, nu, a, bt)
subroutine, public tnsr3d_el_sx(v, nv, u, nu, a, bt, ct)
Tensor operations libxsmm backend.
subroutine, public tnsr3d_el_xsmm(v, nv, u, nu, a, bt, ct)
subroutine, public tnsr3d_xsmm(v, nv, u, nu, a, bt, ct, nelv)
subroutine, public tnsr2d_el_xsmm(v, nv, u, nu, a, bt)
subroutine, public tnsr1_3d_xsmm(v, nv, nu, a, bt, ct, nelv)
subroutine, public tnsr3d_el(v, nv, u, nu, a, bt, ct)
Tensor product performed on a single element.
subroutine, public trsp1(a, n)
In-place transpose of a square tensor.
subroutine, public addtnsr(s, h1, h2, h3, nx, ny, nz)
Maps and adds to S a tensor product form of the three functions H1,H2,H3. This is a single element ro...
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 trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
subroutine triple_tensor_product_scalar(v, u, nu, hr, hs, ht)
Computes the tensor product . This operation is usually performed for spectral interpolation of a sca...
subroutine, public tnsr3d(v, nv, u, nu, a, bt, ct, nelv)
Tensor product performed on nelv elements.
subroutine, public tnsr1_3d(v, nv, nu, a, bt, ct, nelv)
In place tensor product .
subroutine, public tensr3(v, nv, u, nu, a, bt, ct, w)
Tensor product .
subroutine triple_tensor_product_vector(v, u1, u2, u3, nu, hr, hs, ht)
Computes the tensor product on a vector field . This operation is usually performed for spectral inte...
subroutine, public tnsr2d_el(v, nv, u, nu, a, bt)
Computes .