| 
    Neko 0.9.1
    
   A portable framework for high-order spectral element flow simulations 
   | 
 
Fast diagonalization methods from NEKTON.
Functions/Subroutines | |
| subroutine, public | fd_weights_full (xi, x, n, m, c) | 
| Compute finite-difference stencil weights for evaluating derivatives up to order \(m\) at a point.   | |
| 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 = convection operator b*d d = derivative matrix dgll = derivative matrix, mapping from pressure nodes to velocity jgll = interpolation matrix, mapping from pressure nodes to velocity z = GLL points.   | |
| 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 fast3d::fd_weights_full | ( | real(kind=rp), intent(in) | xi, | 
| real(kind=rp), dimension(0:n), intent(in) | x, | ||
| integer, intent(in) | n, | ||
| integer, intent(in) | m, | ||
| real(kind=rp), dimension(0:n,0:m), intent(out) | c | ||
| ) | 
This routine comes from the Appendix C of "A Practical Guide to Pseudospectral Methods" by B. Fornberg, Cambridge University Press, 1996.
Given gridpoints \( x_0, x_1, \dots x_n \) and some point \(\xi\) (not necessarily a grid point!) find weights \( c_{j, k} \), such that the expansions \( \frac{d^k f}{d x^k}|_{x=\xi} \approx \sum_{j=0}^n c_{j,k} f(x_j)\), \(k=0, \dots m\) are optimal. Note that finite-difference stencils are exactly such type of expansions. For the derivation of the algorithm, refer to 3.1 in the reference above.
_full refers to the fact that we use the values \(f(x_j)\) at all available nodes \(x\) to construct the expansion. So we always get the finite difference stencil of maximum order possible.| xi | Point at which the approximations are to be accurate | |
| x | The coordinates for the grid points | |
| [in] | n | The size of x is n + 1  | 
| [in] | m | Highest order of derivative to be approximated | 
| c | The stencil weights. Row j corresponds to weight of \(f(x_j)\) and column k to the kth derivative | 
Definition at line 104 of file fast3d.f90.

| subroutine, public fast3d::semhat | ( | real(kind=rp), dimension(0:n,0:n), intent(inout) | a, | 
| real(kind=rp), dimension(0:n), intent(inout) | b, | ||
| real(kind=rp), dimension(0:n,0:n), intent(inout) | c, | ||
| real(kind=rp), dimension(0:n,0:n), intent(inout) | d, | ||
| real(kind=rp), dimension(0:n), intent(inout) | z, | ||
| real(kind=rp), dimension(0:n,1:n-1), intent(inout) | dgll, | ||
| real(kind=rp), dimension(0:n,1:n-1), intent(inout) | jgll, | ||
| real(kind=rp), dimension(1:n-1), intent(inout) | bgl, | ||
| real(kind=rp), dimension(1:n-1), intent(inout) | zgl, | ||
| real(kind=rp), dimension(1:n-1,0:n), intent(inout) | dgl, | ||
| real(kind=rp), dimension(1:n-1,0:n), intent(inout) | jgl, | ||
| integer, intent(in) | n, | ||
| real(kind=rp), dimension(0:2*n+1), intent(inout) | w | ||
| ) | 
zgl = GL points bgl = diagonal mass matrix on GL dgl = derivative matrix, mapping from velocity nodes to pressure jgl = interpolation matrix, mapping from velocity nodes to pressure
n = polynomial degree (velocity space) w = work array of size 2*n+2
Currently, this is set up for pressure nodes on the interior GLL pts.
Definition at line 167 of file fast3d.f90.


| subroutine, public fast3d::setup_intp | ( | real(kind=rp), dimension(n_to, n_from), intent(inout) | jh, | 
| real(kind=rp), dimension(n_from, n_to), intent(inout) | jht, | ||
| real(kind=rp), dimension(n_to), intent(inout) | z_to, | ||
| real(kind=rp), dimension(n_from), intent(inout) | z_from, | ||
| integer, intent(in) | n_to, | ||
| integer, intent(in) | n_from, | ||
| integer, intent(in) | derivative | ||
| ) | 
This is essentially a wrapper for calling fd_weights_full() for several points. For each point in z_to, we get a set of interpolation weights of size n_from. The result is thus a matrix of weights, each row corresponding to a point in z_to and each column the weight of a point in z_from.
This routine is used for interpolating between elements of different polynomial order. In other words, belonging to different space::space_t . The points are then GL, GLL, etc., depending on the space.
| jh | Matrix of the interpolation weights. | 
| jht | Same as jh but transposed.  | 
| z_to | Target points for interpolation. | 
| z_from | Quadrature points. | 
| n_to | Number of points in z_to.  | 
| n_from | Number of points in z_from.  | 
| derivative | Specifies if we want the derivative interpolation instead, e.g. derivative = 1 refers to the first derivative etc.  | 
Definition at line 242 of file fast3d.f90.

