Line data Source code
1 0 : // Distributed under the MIT License. 2 : // See LICENSE.txt for details. 3 : 4 : #pragma once 5 : 6 : #include <cstddef> 7 : 8 : /// \cond 9 : class Matrix; 10 : template <size_t> 11 : class Mesh; 12 : namespace Spectral { 13 : enum class Basis : uint8_t; 14 : enum class Quadrature : uint8_t; 15 : enum class Parity : uint8_t; 16 : } // namespace Spectral 17 : /// \endcond 18 : 19 : namespace Spectral { 20 : /// @{ 21 : /*! 22 : * \brief %Matrix used to compute the derivative of a function. 23 : * 24 : * \details For a function represented by the nodal coefficients \f$u_j\f$ a 25 : * matrix multiplication with the differentiation matrix \f$D_{ij}\f$ gives the 26 : * coefficients of the function's derivative. Since \f$u(x)\f$ is expanded in 27 : * Lagrange polynomials \f$u(x)=\sum_j u_j l_j(x)\f$ the differentiation matrix 28 : * is computed as \f$D_{ij}=l_j^\prime(\xi_i)\f$ where the \f$\xi_i\f$ are the 29 : * collocation points. 30 : * 31 : * The finite difference matrix uses summation by parts operators, 32 : * \f$D_{2-1}, D_{4-2}, D_{4-3}\f$, and \f$D_{6-5}\f$ from \cite Diener2005tn. 33 : * 34 : * \param num_points The number of collocation points 35 : */ 36 : template <Basis BasisType, Quadrature QuadratureType> 37 1 : const Matrix& differentiation_matrix(size_t num_points); 38 : template <Basis BasisType, Quadrature QuadratureType> 39 1 : const Matrix& differentiation_matrix_transpose(size_t num_points); 40 : /// @} 41 : 42 : /*! 43 : * \brief %Matrix used to compute the derivative of a function of known parity 44 : * 45 : * \details For the Zernike basis, defined on \f$[0,1]\f$, `GaussRadauUpper` 46 : * quadratrue shifts the collocation points to the upper side, which 47 : * contributes to inaccurate differentiation at the lower side due to the low 48 : * density of points. By knowing the parity of functions in this basis as it 49 : * has two indices, we can extend the function to the negative \f$r\f$, 50 : * greatly reducing errors. 51 : * 52 : * \param num_points The number of collocation points 53 : * \param parity The Parity of the function 54 : */ 55 : template <Basis BasisType, Quadrature QuadratureType> 56 1 : const Matrix& differentiation_matrix(size_t num_points, Parity parity); 57 : 58 : /// @{ 59 : /*! 60 : * \brief Differentiation matrix for a one-dimensional mesh. 61 : * 62 : * \see differentiation_matrix(size_t) 63 : */ 64 1 : const Matrix& differentiation_matrix(const Mesh<1>& mesh); 65 1 : const Matrix& differentiation_matrix_transpose(const Mesh<1>& mesh); 66 : /// @} 67 : 68 : /*! 69 : * \brief %Matrix used to compute the divergence of the flux in weak form. 70 : * 71 : * This is the transpose of the differentiation matrix multiplied by quadrature 72 : * weights that appear in DG integrals: 73 : * 74 : * \begin{equation} 75 : * \frac{D^T_{ij}} \frac{w_j}{w_i} 76 : * \end{equation} 77 : * 78 : * \param num_points The number of collocation points 79 : */ 80 : template <Basis BasisType, Quadrature QuadratureType> 81 1 : const Matrix& weak_flux_differentiation_matrix(size_t num_points); 82 : 83 : /*! 84 : * \brief %Matrix used to compute the divergence of the flux in weak form. 85 : * 86 : * \see weak_flux_differentiation_matrix(size_t) 87 : */ 88 1 : const Matrix& weak_flux_differentiation_matrix(const Mesh<1>& mesh); 89 : } // namespace Spectral