Doxygen API reference documentation for ideal.II
|
A class for dG elements in time and arbitrary elements in space. More...
#include <fe_dg.hh>
Public Types | |
enum | support_type { Legendre , Lobatto , RadauLeft , RadauRight } |
Choice of underlying temporal support points. More... | |
Public Member Functions | |
DG_FiniteElement (std::shared_ptr< dealii::FiniteElement< dim >> fe_space, const unsigned int r, support_type type=support_type::Lobatto) | |
Constructor for the finite element class. | |
std::shared_ptr< dealii::FiniteElement< dim > > | spatial () |
The underlying spatial finite element. More... | |
std::shared_ptr< dealii::FiniteElement< 1 > > | temporal () |
The underlying temporal finite element. More... | |
support_type | type () |
The support_type used for construction. More... | |
Public Attributes | |
const unsigned int | dofs_per_cell |
The number of degrees of freedom per space-time element. More... | |
A class for dG elements in time and arbitrary elements in space.
This class implements a tensor product finite element of a discontinuous Galerkin temporal element and a user supplied spatial element. For the temporal element both GaussLegendre and GaussLobatto support points can be chosen
enum idealii::spacetime::DG_FiniteElement::support_type |
Choice of underlying temporal support points.
Allows the choice between GaussLobatto and GaussLegendre support points for the dG(r) elements in time.
The choice of Lobatto results in 1D dG elements as defined in FE_DGQ.
The choice of Legendre leads to dG elements without support points on the interval/element edges. This might lead to better convergence, but the resulting off-diagonal sparsity pattern will be larger as jump terms then depend on all temporal dofs.
For r=0 the resulting element always uses the interval midpoint.
Enumerator | |
---|---|
Legendre | Support points based on QGauss<1> |
Lobatto | for dG(r), r>0: Support points based on QGaussLobatto<1> |
RadauLeft | Support points based on left QGaussRadau<1> |
RadauRight | Support points based on right QGaussRadau<1> |
std::shared_ptr<dealii::FiniteElement<dim> > idealii::spacetime::DG_FiniteElement< dim >::spatial | ( | ) |
The underlying spatial finite element.
std::shared_ptr<dealii::FiniteElement<1> > idealii::spacetime::DG_FiniteElement< dim >::temporal | ( | ) |
The underlying temporal finite element.
support_type idealii::spacetime::DG_FiniteElement< dim >::type | ( | ) |
The support_type used for construction.
This function is mostly used internally. One example would be the construction of the upstream and downstream sparsity patterns, which include off diagonals for jump terms.
const unsigned int idealii::spacetime::DG_FiniteElement< dim >::dofs_per_cell |
The number of degrees of freedom per space-time element.
In practice this is simply (r+1) times the number of DoFs per spatial element.