Represents a basis of Lagrange polynomials. More...
#include <Basis.hpp>
Public Member Functions | |
| Basis (int row_size_arg) | |
| virtual double | node (int i) const =0 |
the ith interpolation node (i.e. quadrature point). Should be in interval [0, 1]. | |
| Mat | nodes () const |
| Get vector of nodes. | |
| virtual Mat | node_weights () const =0 |
node_weights()(i) is the quadrature weight associated with node(i) Integrates in domain [0, 1]. | |
| virtual Mat< dyn, dyn > | diff_mat () const =0 |
| Differentiation matrix. | |
| virtual Mat< dyn, dyn > | boundary () const =0 |
boundary()(i, j) is the jth basis polynomial evaluated at i for i in {0, 1}. | |
| virtual Mat | orthogonal (int degree) const =0 |
orthogonal(degree)(i) is the degreeth-degree Legendre polynomial evaluated at node(i). | |
| virtual Mat< dyn, dyn > | filter () const =0 |
| Matrix that applies a low-pass filter to the Legendre modes of a polynomial. | |
| Mat< dyn, dyn > | interpolate (const Mat<> &sample) const |
| Interpolate to arbitrary points. | |
| virtual Mat< dyn, dyn > | prolong (int i_half) const |
| Matrix to prolong into polynomial space of one higher refinement level. | |
| virtual Mat< dyn, dyn > | restrict (int i_half) const =0 |
| Restrict from refined space above to polynomial space spanned by this basis. | |
| double | max_cfl () const |
| maximum stable CFL number for 1D convection | |
| double | step_ratio () const |
| ratio of time step for convection stage 2 to stage 1 | |
| virtual double | min_eig_diffusion () const =0 |
| nondimensional minimum eigenvalue of 1D diffusion operator | |
Public Attributes | |
| const int | row_size |
| row size | |
Protected Member Functions | |
| virtual double | min_eig_convection () const =0 |
| nondimensional minimum real part of eigenvals of 1D convection operator | |
| virtual double | quadratic_safety () const =0 |
Represents a basis of Lagrange polynomials.
Interpolates between a set of quadrature points, a.k.a. nodes. The ith basis vector is the polynomial with a value of 1 at the ith node and 0 at all other nodes.
| hexed::Basis::Basis | ( | int | row_size_arg | ) |
| row_size_arg | specifies member row_size. |
boundary()(i, j) is the jth basis polynomial evaluated at i for i in {0, 1}.
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.
Differentiation matrix.
diff_mat()(i, j) is the derivative of interpolating polynomial j evaluated at node(i).
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.
Matrix that applies a low-pass filter to the Legendre modes of a polynomial.
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.
Interpolate to arbitrary points.
Current implementation is not very performance-optimized. interpolate(sample)(i, j) is the jth basis polynomial evaluated at sample(i)
|
protectedpure virtual |
nondimensional minimum real part of eigenvals of 1D convection operator
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.
|
pure virtual |
nondimensional minimum eigenvalue of 1D diffusion operator
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.
|
pure virtual |
the ith interpolation node (i.e. quadrature point). Should be in interval [0, 1].
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.
|
pure virtual |
node_weights()(i) is the quadrature weight associated with node(i) Integrates in domain [0, 1].
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.
|
pure virtual |
orthogonal(degree)(i) is the degreeth-degree Legendre polynomial evaluated at node(i).
Polynomial is scaled so that its norm is 1 with respect to the quadrature rule of the basis.
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.
Matrix to prolong into polynomial space of one higher refinement level.
prolong(i_half)(i, j) is the jth basis polynomial evaluated at the ith node in the refined space. If i_half is 0, refined space is interval [0, 0.5]. If i_half is 1, then [0.5, 1]. This is a basic, naive implementation. Derived classes may override for performance and/or precision.
Reimplemented in hexed::Gauss_legendre, and hexed::Gauss_lobatto.
|
protectedpure virtual |
Implemented in hexed::Equidistant.
Restrict from refined space above to polynomial space spanned by this basis.
restrict(i_half)(i, j) is the value at the ith standard node of the orthogonal projection of the jth polynomial in the refined space into the space of this basis.
Implemented in hexed::Equidistant, hexed::Gauss_legendre, and hexed::Gauss_lobatto.