105    integer, 
intent(in) :: n
 
  106    integer, 
intent(in) :: m
 
  107    real(kind=
rp), 
intent(in) :: x(0:n)
 
  108    real(kind=
rp), 
intent(out) :: c(0:n,0:m)
 
  109    real(kind=
rp), 
intent(in) :: xi
 
  110    real(kind=
xp) :: c1, c2, c3, c4, c5
 
  111    integer :: i, j, k, mn
 
  133             c(i,k) = c1 * (k * c(i-1,k-1) - c5 * c(i-1,k)) / c2
 
  135          c(i,0) = -c1 * c5 * c(i-1,0) / c2
 
  137             c(j,k) = (c4 * c(j,k) - k * c(j,k-1)) / c3
 
  139          c(j,0) = c4 * c(j,0) / c3
 
 
  167  subroutine semhat(a, b, c, d, z, dgll, jgll, bgl, zgl, dgl, jgl, n, w)
 
  168    integer, 
intent(in) :: n
 
  169    real(kind=
rp), 
intent(inout) :: a(0:n,0:n)
 
  170    real(kind=
rp), 
intent(inout) :: b(0:n)
 
  171    real(kind=
rp), 
intent(inout) :: c(0:n,0:n)
 
  172    real(kind=
rp), 
intent(inout) :: d(0:n,0:n)
 
  173    real(kind=
rp), 
intent(inout) :: z(0:n)
 
  174    real(kind=
rp), 
intent(inout) :: 
dgll(0:n,1:n-1),jgll(0:n,1:n-1)
 
  175    real(kind=
rp), 
intent(inout) :: bgl(1:n-1)
 
  176    real(kind=
rp), 
intent(inout) :: zgl(1:n-1)
 
  177    real(kind=
rp), 
intent(inout) :: dgl(1:n-1,0:n)
 
  178    real(kind=
rp), 
intent(inout) :: jgl(1:n-1,0:n)
 
  179    real(kind=
rp), 
intent(inout) :: w(0:2*n+1)
 
  180    integer :: np, nm, n2, i, j, k
 
  205             a(i,j) = a(i,j) + d(k,i)*b(k)*d(k,j)
 
  210    call zwgl(zgl, bgl, nm)
 
 
  242  subroutine setup_intp(jh, jht, z_to, z_from, n_to, n_from, derivative)
 
  244    integer, 
intent(in) :: n_to, n_from, 
derivative 
  245    real(kind=
rp), 
intent(inout) :: jh(n_to, n_from), jht(n_from, n_to)
 
  246    real(kind=
rp), 
intent(inout) :: z_to(n_to), z_from(n_from)
 
 
Implements the derivative_t type.
 
Fast diagonalization methods from NEKTON.
 
subroutine, public semhat(a, b, c, d, z, dgll, jgll, bgl, zgl, dgl, jgl, n, w)
Generate matrices for single element, 1D operators: a = Laplacian b = diagonal mass matrix c = convec...
 
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.
 
subroutine, public rzero(a, n)
Zero a real vector.
 
integer, parameter, public xp
 
integer, parameter, public rp
Global precision used in computations.
 
LIBRARY ROUTINES FOR SPECTRAL METHODS.
 
subroutine dgll(d, dt, z, nz, nzd)
 
subroutine zwgll(z, w, np)
 
subroutine zwgl(z, w, np)
Generate NP Gauss Legendre points Z and weights W associated with Jacobi polynomial ....