Neko  0.8.99
A portable framework for high-order spectral element flow simulations
fdm_xsmm.f90
Go to the documentation of this file.
1 
2 module fdm_xsmm
3  use num_types
4  use tensor_xsmm
5  implicit none
6 
7 contains
8 
9  subroutine fdm_do_fast_xsmm(e, r, s, d, nl, ldim, nelv)
10  integer, intent(in) :: nl, nelv, ldim
11  real(kind=rp), intent(inout) :: e(nl**ldim, nelv)
12  real(kind=rp), intent(inout) :: r(nl**ldim, nelv)
13  real(kind=rp), intent(inout) :: s(nl*nl,2,ldim, nelv)
14  real(kind=rp), intent(inout) :: d(nl**ldim, nelv)
15  integer :: ie, nn, i
16 
17  nn = nl**ldim
18  if(.not. ldim .eq. 3) then
19  do ie = 1, nelv
20  call tnsr2d_el_xsmm(e(1,ie), nl, r(1,ie), nl, s(1,2,1,ie), s(1,1,2,ie))
21  do i = 1, nn
22  r(i,ie) = d(i,ie) * e(i,ie)
23  end do
24  call tnsr2d_el_xsmm(e(1,ie), nl, r(1,ie), nl, s(1,1,1,ie), s(1,2,2,ie))
25  end do
26  else
27  do ie = 1, nelv
28  call tnsr3d_el_xsmm(e(1,ie), nl, r(1,ie), nl, &
29  s(1,2,1,ie), s(1,1,2,ie), s(1,1,3,ie))
30  do i = 1, nn
31  r(i,ie) = d(i,ie) * e(i,ie)
32  end do
33  call tnsr3d_el_xsmm(e(1,ie), nl, r(1,ie), nl, &
34  s(1,1,1,ie), s(1,2,2,ie), s(1,2,3,ie))
35  end do
36  end if
37  end subroutine fdm_do_fast_xsmm
38 
39 end module fdm_xsmm
Fast Diagonalization libxsmm backend.
Definition: fdm_xsmm.f90:2
subroutine fdm_do_fast_xsmm(e, r, s, d, nl, ldim, nelv)
Definition: fdm_xsmm.f90:10
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Tensor operations libxsmm backend.
Definition: tensor_xsmm.F90:61
subroutine, public tnsr3d_el_xsmm(v, nv, u, nu, A, Bt, Ct)
Definition: tensor_xsmm.F90:83
subroutine, public tnsr2d_el_xsmm(v, nv, u, nu, A, Bt)
Definition: tensor_xsmm.F90:72