DUBeat  1.0.1
High-order discontinuous Galerkin methods and applications to cardiac electrophysiology
VolumeHandlerDG< dim > Class Template Reference

Class for the main operations on a discontinuous Galerkin volume element. More...

#include <volume_handler_DG.hpp>

Inheritance diagram for VolumeHandlerDG< dim >:

Public Member Functions

 VolumeHandlerDG (const unsigned int degree)
 Constructor. More...
 
 VolumeHandlerDG (VolumeHandlerDG< dim > &VolumeHandlerDG)=default
 Default copy constructor. More...
 
 VolumeHandlerDG (const VolumeHandlerDG< dim > &VolumeHandlerDG)=default
 Default const copy constructor. More...
 
 VolumeHandlerDG (VolumeHandlerDG< dim > &&VolumeHandlerDG)=default
 Default move constructor. More...
 
void reinit (const typename dealii::DoFHandler< dim >::active_cell_iterator &new_cell)
 Reinit objects on the current new_cell. More...
 
virtual dealii::Point< dim > quadrature_real (const unsigned int q) const
 Return the \(q\)-th spatial quadrature point position on the actual element. More...
 
virtual dealii::Point< dim > quadrature_ref (const unsigned int q) const
 Return the \(q\)-th spatial quadrature point position on the reference element. More...
 
virtual double quadrature_weight (const unsigned int q) const
 Return the quadrature weight associated to the \(q\)-th quadrature point. More...
 
dealii::Tensor< 2, dim > get_jacobian_inverse () const
 Inverse of the Jacobian of the reference-to-actual transformation. More...
 
unsigned int get_n_quad_points () const
 Get the number of quadrature points on the current element. More...
 
virtual ~VolumeHandlerDG ()=default
 Destructor. More...
 

Protected Attributes

const unsigned int n_quad_points_1D
 Number of quadrature points in one dimensional elements. More...
 
dealii::DoFHandler< dim >::active_cell_iterator cell
 Actual DG cell. More...
 
const std::unique_ptr< dealii::FE_SimplexDGP< dim > > fe_dg
 Internal Lagrangian basis class. More...
 
const std::unique_ptr< dealii::MappingFE< dim > > mapping
 Mapping of the discretized space, needed for geometrical reference-to-actual operations. More...
 
const QGaussLegendreSimplex< dim > QGLpoints
 Quadrature formula for volume elements. More...
 
std::unique_ptr< dealii::FEValues< dim > > fe_values
 Internal FEM basis class. More...
 
bool initialized = false
 A condition to inform if the class is initialized on an element or not. More...
 

Detailed Description

template<unsigned int dim>
class VolumeHandlerDG< dim >

Class for the main operations on a discontinuous Galerkin volume element.

Definition at line 50 of file volume_handler_DG.hpp.

Constructor & Destructor Documentation

◆ VolumeHandlerDG() [1/4]

template<unsigned int dim>
VolumeHandlerDG< dim >::VolumeHandlerDG ( const unsigned int  degree)
inline

Constructor.

Definition at line 80 of file volume_handler_DG.hpp.

◆ VolumeHandlerDG() [2/4]

template<unsigned int dim>
VolumeHandlerDG< dim >::VolumeHandlerDG ( VolumeHandlerDG< dim > &  VolumeHandlerDG)
default

Default copy constructor.

◆ VolumeHandlerDG() [3/4]

template<unsigned int dim>
VolumeHandlerDG< dim >::VolumeHandlerDG ( const VolumeHandlerDG< dim > &  VolumeHandlerDG)
default

Default const copy constructor.

◆ VolumeHandlerDG() [4/4]

template<unsigned int dim>
VolumeHandlerDG< dim >::VolumeHandlerDG ( VolumeHandlerDG< dim > &&  VolumeHandlerDG)
default

Default move constructor.

◆ ~VolumeHandlerDG()

template<unsigned int dim>
virtual VolumeHandlerDG< dim >::~VolumeHandlerDG ( )
virtualdefault

Destructor.

Member Function Documentation

◆ get_jacobian_inverse()

template<unsigned int dim>
dealii::Tensor< 2, dim > VolumeHandlerDG< dim >::get_jacobian_inverse

Inverse of the Jacobian of the reference-to-actual transformation.

Definition at line 189 of file volume_handler_DG.hpp.

◆ get_n_quad_points()

template<unsigned int dim>
unsigned int VolumeHandlerDG< dim >::get_n_quad_points

Get the number of quadrature points on the current element.

Definition at line 212 of file volume_handler_DG.hpp.

◆ quadrature_real()

template<unsigned int dim>
dealii::Point< dim > VolumeHandlerDG< dim >::quadrature_real ( const unsigned int  q) const
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.

◆ quadrature_ref()

template<unsigned int dim>
dealii::Point< dim > VolumeHandlerDG< dim >::quadrature_ref ( const unsigned int  q) const
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.

◆ quadrature_weight()

template<unsigned int dim>
double VolumeHandlerDG< dim >::quadrature_weight ( const unsigned int  q) const
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.

◆ reinit()

template<unsigned int dim>
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.

Member Data Documentation

◆ cell

template<unsigned int dim>
dealii::DoFHandler<dim>::active_cell_iterator VolumeHandlerDG< dim >::cell
protected

Actual DG cell.

Definition at line 58 of file volume_handler_DG.hpp.

◆ fe_dg

template<unsigned int dim>
const std::unique_ptr<dealii::FE_SimplexDGP<dim> > VolumeHandlerDG< dim >::fe_dg
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.

◆ fe_values

template<unsigned int dim>
std::unique_ptr<dealii::FEValues<dim> > VolumeHandlerDG< dim >::fe_values
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.

◆ initialized

template<unsigned int dim>
bool VolumeHandlerDG< dim >::initialized = false
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.

◆ mapping

template<unsigned int dim>
const std::unique_ptr<dealii::MappingFE<dim> > VolumeHandlerDG< dim >::mapping
protected

Mapping of the discretized space, needed for geometrical reference-to-actual operations.

Definition at line 68 of file volume_handler_DG.hpp.

◆ n_quad_points_1D

template<unsigned int dim>
const unsigned int VolumeHandlerDG< dim >::n_quad_points_1D
protected

Number of quadrature points in one dimensional elements.

Default is polynomial degree + 2.

Definition at line 55 of file volume_handler_DG.hpp.

◆ QGLpoints

template<unsigned int dim>
const QGaussLegendreSimplex<dim> VolumeHandlerDG< dim >::QGLpoints
protected

Quadrature formula for volume elements.

Definition at line 71 of file volume_handler_DG.hpp.


The documentation for this class was generated from the following file: