Functionality associated with a particular choice of basis functions and quadrature for spectral operations.
More...
|
const Matrix & | projection_matrix_mortar_to_element (MortarSize size, const Mesh< 1 > &element_mesh, const Mesh< 1 > &mortar_mesh) noexcept |
| The projection matrix from a 1D mortar to an element. More...
|
|
const Matrix & | projection_matrix_element_to_mortar (MortarSize size, const Mesh< 1 > &mortar_mesh, const Mesh< 1 > &element_mesh) noexcept |
| The projection matrix from a 1D element to a mortar. See projection_matrix_mortar_to_element() for details.
|
|
std::ostream & | operator<< (std::ostream &os, const Basis &basis) noexcept |
|
std::ostream & | operator<< (std::ostream &os, const Quadrature &quadrature) noexcept |
|
template<Basis BasisType, Quadrature QuadratureType> |
const DataVector & | collocation_points (size_t num_points) noexcept |
| Collocation points. More...
|
|
template<Basis BasisType, Quadrature QuadratureType> |
const DataVector & | quadrature_weights (size_t num_points) noexcept |
| Weights to compute definite integrals. More...
|
|
template<Basis BasisType, Quadrature QuadratureType, typename T > |
Matrix | interpolation_matrix (size_t num_points, const T &target_points) noexcept |
| Matrix used to interpolate to the target_points . More...
|
|
template<typename T > |
Matrix | interpolation_matrix (const Mesh< 1 > &mesh, const T &target_points) noexcept |
| Interpolation matrix to the target_points for a one-dimensional mesh. More...
|
|
template<Basis BasisType, typename T > |
T | compute_basis_function_value (size_t k, const T &x) noexcept |
| Compute the function values of the basis function \(\Phi_k(x)\) (zero-indexed).
|
|
template<Basis > |
DataVector | compute_inverse_weight_function_values (const DataVector &) noexcept |
| Compute the inverse of the weight function \(w(x)\) w.r.t. which the basis functions are orthogonal. See the description of quadrature_weights(size_t) for details.
|
|
template<Basis BasisType> |
double | compute_basis_function_normalization_square (size_t k) noexcept |
| Compute the normalization square of the basis function \(\Phi_k\) (zero-indexed), i.e. the weighted definite integral over its square.
|
|
template<Basis BasisType, Quadrature QuadratureType> |
std::pair< DataVector, DataVector > | compute_collocation_points_and_weights (size_t num_points) noexcept |
| Compute the collocation points and weights associated to the basis and quadrature. More...
|
|
const DataVector & | collocation_points (const Mesh< 1 > &mesh) noexcept |
| Collocation points for a one-dimensional mesh. More...
|
|
const DataVector & | quadrature_weights (const Mesh< 1 > &mesh) noexcept |
| Quadrature weights for a one-dimensional mesh. More...
|
|
template<Basis BasisType, Quadrature QuadratureType> |
const Matrix & | differentiation_matrix (size_t num_points) noexcept |
| Matrix used to compute the derivative of a function. More...
|
|
const Matrix & | differentiation_matrix (const Mesh< 1 > &mesh) noexcept |
| Differentiation matrix for a one-dimensional mesh. More...
|
|
template<Basis BasisType, Quadrature QuadratureType> |
const Matrix & | spectral_to_grid_points_matrix (size_t num_points) noexcept |
| Matrix used to transform from the spectral coefficients (modes) of a function to its nodal coefficients. Also referred to as the Vandermonde matrix. More...
|
|
const Matrix & | spectral_to_grid_points_matrix (const Mesh< 1 > &mesh) noexcept |
| Transformation matrix from modal to nodal coefficients for a one-dimensional mesh. More...
|
|
template<Basis BasisType, Quadrature QuadratureType> |
const Matrix & | grid_points_to_spectral_matrix (size_t num_points) noexcept |
| Matrix used to transform from the nodal coefficients of a function to its spectral coefficients (modes). Also referred to as the inverse Vandermonde matrix. More...
|
|
const Matrix & | grid_points_to_spectral_matrix (const Mesh< 1 > &mesh) noexcept |
| Transformation matrix from nodal to modal coefficients for a one-dimensional mesh. More...
|
|
template<Basis BasisType, Quadrature QuadratureType> |
const Matrix & | linear_filter_matrix (size_t num_points) noexcept |
| Matrix used to linearize a function. More...
|
|
const Matrix & | linear_filter_matrix (const Mesh< 1 > &mesh) noexcept |
| Linear filter matrix for a one-dimensional mesh. More...
|
|
Functionality associated with a particular choice of basis functions and quadrature for spectral operations.
Details
The functions in this namespace provide low-level access to collocation points, quadrature weights and associated matrices, such as for differentiation and interpolation. They are available in two versions: either templated directly on the enum cases of the Spectral::Basis and Spectral::Quadrature types, or taking a one-dimensional Mesh as their argument.
- Note
- Generally you should prefer working with a Mesh. Use its Mesh::slice_through() method to retrieve the mesh in a particular dimension:
{{3, 4}},
{{Spectral::Basis::Legendre, Spectral::Basis::Legendre}},
{{Spectral::Quadrature::Gauss, Spectral::Quadrature::GaussLobatto}}};
const auto collocation_points_in_first_dim =
Most algorithms in this namespace are adapted from [16].
template<Basis BasisType, Quadrature QuadratureType>
const DataVector & Spectral::quadrature_weights |
( |
size_t |
num_points | ) |
|
|
noexcept |
Weights to compute definite integrals.
Details
These are the coefficients to contract with the nodal function values \(f_k\) to approximate the definite integral \(I[f]=\int f(x)\mathrm{d}x\).
Note that the term quadrature also often refers to the quantity \(Q[f]=\int f(x)w(x)\mathrm{d}x\approx \sum_k f_k w_k\). Here, \(w(x)\) denotes the basis-specific weight function w.r.t. to which the basis functions \(\Phi_k\) are orthogonal, i.e \(\int\Phi_i(x)\Phi_j(x)w(x)=0\) for \(i\neq j\). The weights \(w_k\) approximate this inner product. To approximate the definite integral \(I[f]\) we must employ the coefficients \(\frac{w_k}{w(\xi_k)}\) instead, where the \(\xi_k\) are the collocation points. These are the coefficients this function returns. Only for a unit weight function \(w(x)=1\), i.e. a Legendre basis, is \(I[f]=Q[f]\) so this function returns the \(w_k\) identically.
- Parameters
-
num_points | The number of collocation points |