|
DUBeat 1.0.1
High-order discontinuous Galerkin methods and applications to cardiac electrophysiology
|
Class for the main operations on a discontinuous Galerkin volume element. More...
#include <volume_handler_DG.hpp>

Public Member Functions | |
| VolumeHandlerDG (const unsigned int degree) | |
| Constructor. | |
| VolumeHandlerDG (VolumeHandlerDG< dim > &VolumeHandlerDG)=default | |
| Default copy constructor. | |
| VolumeHandlerDG (const VolumeHandlerDG< dim > &VolumeHandlerDG)=default | |
| Default const copy constructor. | |
| VolumeHandlerDG (VolumeHandlerDG< dim > &&VolumeHandlerDG)=default | |
| Default move constructor. | |
| void | reinit (const typename dealii::DoFHandler< dim >::active_cell_iterator &new_cell) |
| Reinit objects on the current new_cell. | |
| virtual dealii::Point< dim > | quadrature_real (const unsigned int q) const |
| Return the \(q\)-th spatial quadrature point position on the actual element. | |
| virtual dealii::Point< dim > | quadrature_ref (const unsigned int q) const |
| Return the \(q\)-th spatial quadrature point position on the reference element. | |
| virtual double | quadrature_weight (const unsigned int q) const |
| Return the quadrature weight associated to the \(q\)-th quadrature point. | |
| dealii::Tensor< 2, dim > | get_jacobian_inverse () const |
| Inverse of the Jacobian of the reference-to-actual transformation. | |
| unsigned int | get_n_quad_points () const |
| Get the number of quadrature points on the current element. | |
| virtual | ~VolumeHandlerDG ()=default |
| Destructor. | |
Protected Attributes | |
| const unsigned int | n_quad_points_1D |
| Number of quadrature points in one dimensional elements. | |
| dealii::DoFHandler< dim >::active_cell_iterator | cell |
| Actual DG cell. | |
| const std::unique_ptr< dealii::FE_SimplexDGP< dim > > | fe_dg |
| Internal Lagrangian basis class. | |
| const std::unique_ptr< dealii::MappingFE< dim > > | mapping |
| Mapping of the discretized space, needed for geometrical reference-to-actual operations. | |
| const QGaussLegendreSimplex< dim > | QGLpoints |
| Quadrature formula for volume elements. | |
| std::unique_ptr< dealii::FEValues< dim > > | fe_values |
| Internal FEM basis class. | |
| bool | initialized = false |
| A condition to inform if the class is initialized on an element or not. | |
Class for the main operations on a discontinuous Galerkin volume element.
Definition at line 50 of file volume_handler_DG.hpp.
|
inline |
Constructor.
Definition at line 80 of file volume_handler_DG.hpp.
|
default |
Default copy constructor.
|
default |
Default const copy constructor.
|
default |
Default move constructor.
|
virtualdefault |
Destructor.
| dealii::Tensor< 2, dim > VolumeHandlerDG< dim >::get_jacobian_inverse | ( | ) | const |
Inverse of the Jacobian of the reference-to-actual transformation.
Definition at line 189 of file volume_handler_DG.hpp.
| unsigned int VolumeHandlerDG< dim >::get_n_quad_points | ( | ) | const |
Get the number of quadrature points on the current element.
Definition at line 212 of file volume_handler_DG.hpp.
|
virtual |
Return the \(q\)-th spatial quadrature point position on the actual element.
Reimplemented in FaceHandlerDG< dim >, and FaceHandlerDG< lifex::dim >.
Definition at line 150 of file volume_handler_DG.hpp.
|
virtual |
Return the \(q\)-th spatial quadrature point position on the reference element.
Reimplemented in FaceHandlerDG< dim >, and FaceHandlerDG< lifex::dim >.
Definition at line 163 of file volume_handler_DG.hpp.
|
virtual |
Return the quadrature weight associated to the \(q\)-th quadrature point.
Reimplemented in FaceHandlerDG< dim >, and FaceHandlerDG< lifex::dim >.
Definition at line 176 of file volume_handler_DG.hpp.
| void VolumeHandlerDG< dim >::reinit | ( | const typename dealii::DoFHandler< dim >::active_cell_iterator & | new_cell | ) |
Reinit objects on the current new_cell.
Definition at line 140 of file volume_handler_DG.hpp.
|
protected |
Actual DG cell.
Definition at line 58 of file volume_handler_DG.hpp.
|
protected |
Internal Lagrangian basis class.
This internal member permits to exploit useful already implemented operations. Polynomial order is always set to 1 because the class executes only geometric operations not related to the degrees of freedom.
Definition at line 64 of file volume_handler_DG.hpp.
|
protected |
Internal FEM basis class.
This internal member permits to exploit useful already implemented operations. As for fe_dg, the polynomial order is always 1 because it is only needed for geometric operations not related to the degrees of freedom.
Definition at line 77 of file volume_handler_DG.hpp.
|
protected |
A condition to inform if the class is initialized on an element or not.
Definition at line 80 of file volume_handler_DG.hpp.
|
protected |
Mapping of the discretized space, needed for geometrical reference-to-actual operations.
Definition at line 68 of file volume_handler_DG.hpp.
|
protected |
Number of quadrature points in one dimensional elements.
Default is polynomial degree + 2.
Definition at line 55 of file volume_handler_DG.hpp.
|
protected |
Quadrature formula for volume elements.
Definition at line 71 of file volume_handler_DG.hpp.