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

Class for the assembly of the main local matrices for discontinuous Galerkin methods. More...

#include <assemble_DG.hpp>

Collaboration diagram for AssembleDG< basis >:

Public Member Functions

 AssembleDG (const unsigned int degree)
 Constructor. More...
 
 AssembleDG (AssembleDG< basis > &)=default
 Default copy constructor. More...
 
 AssembleDG (const AssembleDG< basis > &)=default
 Default const copy constructor. More...
 
 AssembleDG (AssembleDG< basis > &&)=default
 Default move constructor. More...
 
void reinit (const typename dealii::DoFHandler< lifex::dim >::active_cell_iterator &new_cell, const unsigned int new_edge)
 Reinitialize object on the current new_edge of the new_cell. More...
 
void reinit (const typename dealii::DoFHandler< lifex::dim >::active_cell_iterator &new_cell)
 Reinitialize object on the current new_cell. More...
 
unsigned int get_dofs_per_cell () const
 Return the number of degrees of freedom per element. More...
 
dealii::FullMatrix< double > local_V () const
 Assembly of stiffness matrix: More...
 
dealii::FullMatrix< double > local_V (const dealii::Tensor< 2, lifex::dim > &Sigma) const
 Assembly of stiffness matrix with diffusion parameter: More...
 
dealii::FullMatrix< double > local_M () const
 Assembly of local mass matrix: More...
 
dealii::Vector< double > local_rhs (const std::shared_ptr< dealii::Function< lifex::dim >> &f_ex) const
 Assembly of the local forcing term: More...
 
dealii::Vector< double > local_rhs_edge_dirichlet (const double stability_coefficient, const double theta, const std::shared_ptr< lifex::utils::FunctionDirichlet > &u_ex) const
 Assembly of the local matrix associated to non-homogeneous Dirichlet boundary conditions: More...
 
dealii::Vector< double > local_rhs_edge_neumann (const std::shared_ptr< dealii::Function< lifex::dim >> &g) const
 Assembly of the right hand side term associated to non-homogeneous Neumann boundary conditions: More...
 
dealii::FullMatrix< double > local_SC (const double stability_coefficient) const
 Assembly of the component of the local matrix S that is evaluated on the same side of the edge: More...
 
dealii::FullMatrix< double > local_SC (const double stability_coefficient, const dealii::Tensor< 2, lifex::dim > &Sigma) const
 Assembly of the component of the local matrix S that is evaluated on the same side of the edge with diffusion parameter: More...
 
dealii::FullMatrix< double > local_SN (const double stability_coefficient) const
 Assembly of the component of the local matrix S that is evaluated on the two sides of the edge: More...
 
dealii::FullMatrix< double > local_SN (const double stability_coefficient, const dealii::Tensor< 2, lifex::dim > &Sigma) const
 Assembly of the component of the local matrix S that is evaluated on the two sides of the edge with diffusion parameter: More...
 
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > local_IC (const double theta) const
 Assembly of the component of the local matrix I that is evaluated on the same side of the edge. More...
 
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > local_IC (const double theta, const dealii::Tensor< 2, lifex::dim > &Sigma) const
 Assembly of the component of the local matrix I that is evaluated on the same side of the edge with diffusion parameter. More...
 
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > local_IB (const double theta) const
 Assembly of the component of the local matrix I that is evaluated on the boundary edges. More...
 
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > local_IB (const double theta, const dealii::Tensor< 2, lifex::dim > &Sigma) const
 Assembly of the component of the local matrix I that is evaluated on the boundary edges with diffusion parameter. More...
 
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > local_IN (const double theta) const
 Assembly of the component of the local matrix I that is evaluated on the two sides of the edge. More...
 
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > local_IN (const double theta, const dealii::Tensor< 2, lifex::dim > &Sigma) const
 Assembly of the component of the local matrix I that is evaluated on the two sides of the edge with diffusion parameter. More...
 
dealii::Vector< double > local_u0_M_rhs (const lifex::LinAlg::MPI::Vector &u0, const std::vector< dealii::types::global_dof_index > &dof_indices) const
 Assembly of the right hand side term associated to the previous time-step solution in time-dependent problems: More...
 
dealii::Vector< double > local_w0_M_rhs (const lifex::LinAlg::MPI::Vector &u0, const std::vector< dealii::types::global_dof_index > &dof_indices) const
 Assembly of the right hand side term associated to the previous time-step gating variable solution (for monodomain problem): More...
 
dealii::FullMatrix< double > local_non_linear_fitzhugh (const lifex::LinAlg::MPI::Vector &u0, const double a, const std::vector< dealii::types::global_dof_index > &dof_indices) const
 Assembly of the non-linear local matrix of the Fitzhugh-Nagumo model: More...
 
virtual ~AssembleDG ()=default
 Destructor. More...
 

Protected Attributes

const std::unique_ptr< basis > basis_ptr
 Basis function class. More...
 
VolumeHandlerDG< lifex::dim > vol_handler
 Volume element handler. More...
 
FaceHandlerDG< lifex::dim > face_handler
 Face element handler. More...
 
FaceHandlerDG< lifex::dim > face_handler_neigh
 Face element handler of the neighbor face. More...
 
unsigned int dofs_per_cell
 Degree of freedom for each cell. More...
 
const unsigned int n_quad_points_1D
 Number of quadrature points in 1D. More...
 
const unsigned int n_quad_points
 Number of quadrature points in the volume: (n_quad_points_1D)^(dim). More...
 
const unsigned int n_quad_points_face
 Number of quadrature points in the face: (n_quad_points_1D)^(dim-1). More...
 
const unsigned int poly_degree
 Polynomial degree. More...
 
dealii::DoFHandler< lifex::dim >::active_cell_iterator cell
 Cell object. More...
 

Detailed Description

template<class basis>
class AssembleDG< basis >

Class for the assembly of the main local matrices for discontinuous Galerkin methods.

Definition of the main terms:

  • \(\mathcal{K}\) as the volume of a generic element.
  • \(\mathcal{F}\) as a generic face of one element.
  • \(\varphi_i\) as the \(i\)-th basis function of the element.
  • The \(+\) superscript to indicate that a function evaluated on the face is the one associated to the same element.
  • The \(-\) superscript to indicate that a function evaluated on the face is the one associated to the neighbor element.
  • \(n^+, n^-\) as the outward/inward normal vector with respect to one element face ( \(n=n^+\) by default).

Definition at line 68 of file assemble_DG.hpp.

Constructor & Destructor Documentation

◆ AssembleDG() [1/4]

template<class basis >
AssembleDG< basis >::AssembleDG ( const unsigned int  degree)
inline

Constructor.

Definition at line 99 of file assemble_DG.hpp.

◆ AssembleDG() [2/4]

template<class basis >
AssembleDG< basis >::AssembleDG ( AssembleDG< basis > &  )
default

Default copy constructor.

◆ AssembleDG() [3/4]

template<class basis >
AssembleDG< basis >::AssembleDG ( const AssembleDG< basis > &  )
default

Default const copy constructor.

◆ AssembleDG() [4/4]

template<class basis >
AssembleDG< basis >::AssembleDG ( AssembleDG< basis > &&  )
default

Default move constructor.

◆ ~AssembleDG()

template<class basis >
virtual AssembleDG< basis >::~AssembleDG ( )
virtualdefault

Destructor.

Member Function Documentation

◆ get_dofs_per_cell()

template<class basis >
unsigned int AssembleDG< basis >::get_dofs_per_cell

Return the number of degrees of freedom per element.

Definition at line 355 of file assemble_DG.hpp.

◆ local_IB() [1/2]

template<class basis >
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > AssembleDG< basis >::local_IB ( const double  theta) const

Assembly of the component of the local matrix I that is evaluated on the boundary edges.

The method returns the two matrices:

\[\begin{aligned} \theta \cdot IB(i,j)=& \: - \theta \int_{\mathcal{F}} \nabla \varphi_i^{+} \cdot n^{+} \varphi_j^{+} \, ds \\ IB^T(i,j)=& - \int_{\mathcal{F}} \nabla \varphi_j^{+} \cdot n^{+} \varphi_i^{+} \, ds \end{aligned} \]

where \(\theta\) is the penalty coefficient.

Definition at line 803 of file assemble_DG.hpp.

◆ local_IB() [2/2]

template<class basis >
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > AssembleDG< basis >::local_IB ( const double  theta,
const dealii::Tensor< 2, lifex::dim > &  Sigma 
) const

Assembly of the component of the local matrix I that is evaluated on the boundary edges with diffusion parameter.

The method returns the two matrices:

\[\begin{aligned} \theta \cdot IB(i,j)=& \: - \theta \int_{\mathcal{F}} \Sigma \nabla \varphi_i^{+} \cdot n^{+} \varphi_j^{+} \, ds \\ IB^T(i,j)=& - \int_{\mathcal{F}} \Sigma \nabla \varphi_j^{+} \cdot n^{+} \varphi_i^{+} \, ds \end{aligned} \]

where \(\theta\) is the penalty coefficient and \(\Sigma\) is the diffusion tensor.

Definition at line 845 of file assemble_DG.hpp.

◆ local_IC() [1/2]

template<class basis >
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > AssembleDG< basis >::local_IC ( const double  theta) const

Assembly of the component of the local matrix I that is evaluated on the same side of the edge.

The method returns the two matrices:

\[ \begin{aligned} \theta \cdot IC(i,j)= & - \frac{\theta}{2} \int_{\mathcal{F}} \nabla \varphi_i^{+} \cdot n^{+} \varphi_j^{+} \, ds \\ IC^T(i,j)= & - \frac{1}{2} \int_{\mathcal{F}} \nabla \varphi_j^{+} \cdot n^{+} \varphi_i^{+} \, ds \end{aligned} \]

where \(\theta\) is the penalty coefficient.

Definition at line 714 of file assemble_DG.hpp.

◆ local_IC() [2/2]

template<class basis >
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > AssembleDG< basis >::local_IC ( const double  theta,
const dealii::Tensor< 2, lifex::dim > &  Sigma 
) const

Assembly of the component of the local matrix I that is evaluated on the same side of the edge with diffusion parameter.

The method returns the two matrices:

\[ \begin{aligned} \theta \cdot IC(i,j)= & - \frac{\theta}{2} \int_{\mathcal{F}} \Sigma \nabla \varphi_i^{+} \cdot n^{+} \varphi_j^{+} \, ds \\ IC^T(i,j)= & - \frac{1}{2} \int_{\mathcal{F}} \Sigma \nabla \varphi_j^{+} \cdot n^{+} \varphi_i^{+} \, ds \end{aligned} \]

where \(\theta\) is the penalty coefficient and \(\Sigma\) is the diffusion tensor.

Definition at line 757 of file assemble_DG.hpp.

◆ local_IN() [1/2]

template<class basis >
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > AssembleDG< basis >::local_IN ( const double  theta) const

Assembly of the component of the local matrix I that is evaluated on the two sides of the edge.

The method returns the two matrices:

\[\begin{aligned} \theta \cdot IN(i,j)=& \frac{\theta}{2} \int_{\mathcal{F}} \nabla \varphi_i^{+} \cdot n^{+} \varphi_j^{-} \, ds \\ IN^T(i,j)=& \frac{1}{2} \int_{\mathcal{F}} \nabla \varphi_j^{+} \cdot n^{+} \varphi_i^{-} \, ds \end{aligned} \]

where \(\theta\) is the penalty coefficient.

Definition at line 890 of file assemble_DG.hpp.

◆ local_IN() [2/2]

template<class basis >
std::pair< dealii::FullMatrix< double >, dealii::FullMatrix< double > > AssembleDG< basis >::local_IN ( const double  theta,
const dealii::Tensor< 2, lifex::dim > &  Sigma 
) const

Assembly of the component of the local matrix I that is evaluated on the two sides of the edge with diffusion parameter.

The method returns the two matrices:

\[\begin{aligned} \theta \cdot IN(i,j)=& \frac{\theta}{2} \int_{\mathcal{F}} \Sigma \nabla \varphi_i^{+} \cdot n^{+} \varphi_j^{-} \, ds \\ IN^T(i,j)=& \frac{1}{2} \int_{\mathcal{F}} \Sigma \nabla \varphi_j^{+} \cdot n^{+} \varphi_i^{-} \, ds \end{aligned} \]

where \(\theta\) is the penalty coefficient and \(\Sigma\) is the diffusion tensor.

Definition at line 936 of file assemble_DG.hpp.

◆ local_M()

template<class basis >
dealii::FullMatrix< double > AssembleDG< basis >::local_M

Assembly of local mass matrix:

\[M(i,j)=\int_{\mathcal{K}} \varphi_j \varphi_i \, dx \]

Definition at line 435 of file assemble_DG.hpp.

◆ local_non_linear_fitzhugh()

template<class basis >
dealii::FullMatrix< double > AssembleDG< basis >::local_non_linear_fitzhugh ( const lifex::LinAlg::MPI::Vector &  u0,
const double  a,
const std::vector< dealii::types::global_dof_index > &  dof_indices 
) const

Assembly of the non-linear local matrix of the Fitzhugh-Nagumo model:

\[C(i,j)=\int_{\mathcal{K}} (u^n-1)(u^n-a) \varphi_j \varphi_i \, dx \]

where \(a\) is parameter of the monodomain equation and \(u^n\) is the solution at a generic previous step \( n\).

Definition at line 1051 of file assemble_DG.hpp.

◆ local_rhs()

template<class basis >
dealii::Vector< double > AssembleDG< basis >::local_rhs ( const std::shared_ptr< dealii::Function< lifex::dim >> &  f_ex) const

Assembly of the local forcing term:

\[F(i)=\int_{\mathcal{K}} f \varphi_i \, dx \]

where \( f \) is the known forcing term of the problem.

Definition at line 461 of file assemble_DG.hpp.

◆ local_rhs_edge_dirichlet()

template<class basis >
dealii::Vector< double > AssembleDG< basis >::local_rhs_edge_dirichlet ( const double  stability_coefficient,
const double  theta,
const std::shared_ptr< lifex::utils::FunctionDirichlet > &  u_ex 
) const

Assembly of the local matrix associated to non-homogeneous Dirichlet boundary conditions:

\[F(i)=\int_{\mathcal{F}} \gamma u \varphi_i^{+} - \theta u \nabla \varphi_i^{+} \cdot n \,ds \]

where \(\gamma\) is the regularity coeficient, \(\theta\) is the penalty coefficient and \(u\) the known solution of the problem.

Definition at line 489 of file assemble_DG.hpp.

◆ local_rhs_edge_neumann()

template<class basis >
dealii::Vector< double > AssembleDG< basis >::local_rhs_edge_neumann ( const std::shared_ptr< dealii::Function< lifex::dim >> &  g) const

Assembly of the right hand side term associated to non-homogeneous Neumann boundary conditions:

\[F(i)=\int_{\mathcal{F}} \varphi_i g \, ds \]

where \(g\) is the gradient of the known solution of the problem.

Definition at line 540 of file assemble_DG.hpp.

◆ local_SC() [1/2]

template<class basis >
dealii::FullMatrix< double > AssembleDG< basis >::local_SC ( const double  stability_coefficient) const

Assembly of the component of the local matrix S that is evaluated on the same side of the edge:

\[SC(i,j)=\int_{\mathcal{F}} \gamma \varphi_j^{+} \varphi_i^{+} \, ds \]

where \( \gamma \) is the regularity coefficient.

Definition at line 569 of file assemble_DG.hpp.

◆ local_SC() [2/2]

template<class basis >
dealii::FullMatrix< double > AssembleDG< basis >::local_SC ( const double  stability_coefficient,
const dealii::Tensor< 2, lifex::dim > &  Sigma 
) const

Assembly of the component of the local matrix S that is evaluated on the same side of the edge with diffusion parameter:

\[SC(i,j)=\int_{\mathcal{F}} \Gamma \varphi_j^{+} \varphi_i^{+} \, ds \]

where \( \Gamma = \gamma(n^T \Sigma n) \), \(\gamma\) is the regularity coefficient and \(\Sigma\) is the diffusion tensor.

Definition at line 601 of file assemble_DG.hpp.

◆ local_SN() [1/2]

template<class basis >
dealii::FullMatrix< double > AssembleDG< basis >::local_SN ( const double  stability_coefficient) const

Assembly of the component of the local matrix S that is evaluated on the two sides of the edge:

\[ SN(i,j)=- \int_{\mathcal{F}} \gamma \varphi_j^{+} \varphi_i^{-} \, ds \]

Definition at line 638 of file assemble_DG.hpp.

◆ local_SN() [2/2]

template<class basis >
dealii::FullMatrix< double > AssembleDG< basis >::local_SN ( const double  stability_coefficient,
const dealii::Tensor< 2, lifex::dim > &  Sigma 
) const

Assembly of the component of the local matrix S that is evaluated on the two sides of the edge with diffusion parameter:

\[ SN(i,j)=- \int_{\mathcal{F}} \Gamma \varphi_j^{+} \varphi_i^{-} \, ds \]

where \( \Gamma = \gamma(n^T \Sigma n) \), \(\gamma\) is the regularity coefficient and \(\Sigma\) is the diffusion tensor.

Definition at line 674 of file assemble_DG.hpp.

◆ local_u0_M_rhs()

template<class basis >
dealii::Vector< double > AssembleDG< basis >::local_u0_M_rhs ( const lifex::LinAlg::MPI::Vector &  u0,
const std::vector< dealii::types::global_dof_index > &  dof_indices 
) const

Assembly of the right hand side term associated to the previous time-step solution in time-dependent problems:

\[F(i)= \int_{\mathcal{K}} u^n \varphi_i \, dx \]

where \(u^n\) is the solution at a generic previous step \( n\).

Definition at line 985 of file assemble_DG.hpp.

◆ local_V() [1/2]

template<class basis >
dealii::FullMatrix< double > AssembleDG< basis >::local_V

Assembly of stiffness matrix:

\[V(i,j)=\int_{\mathcal{K}} \nabla \varphi_j \cdot \nabla \varphi_i \,dx \]

Definition at line 376 of file assemble_DG.hpp.

◆ local_V() [2/2]

template<class basis >
dealii::FullMatrix< double > AssembleDG< basis >::local_V ( const dealii::Tensor< 2, lifex::dim > &  Sigma) const

Assembly of stiffness matrix with diffusion parameter:

\[V(i,j)=\int_{\mathcal{K}} \nabla \varphi_j \cdot \Sigma \nabla \varphi_i \,dx \]

where \(\Sigma\) is the diffusion tensor.

Definition at line 405 of file assemble_DG.hpp.

◆ local_w0_M_rhs()

template<class basis >
dealii::Vector< double > AssembleDG< basis >::local_w0_M_rhs ( const lifex::LinAlg::MPI::Vector &  u0,
const std::vector< dealii::types::global_dof_index > &  dof_indices 
) const

Assembly of the right hand side term associated to the previous time-step gating variable solution (for monodomain problem):

\[F(i)=\int_{\mathcal{K}} w^n \varphi_i \, dx \]

where \(w^n\) is the gating variable solution at a generic previous step \( n\).

Definition at line 1023 of file assemble_DG.hpp.

◆ reinit() [1/2]

template<class basis >
void AssembleDG< basis >::reinit ( const typename dealii::DoFHandler< lifex::dim >::active_cell_iterator new_cell)

Reinitialize object on the current new_cell.

Definition at line 346 of file assemble_DG.hpp.

◆ reinit() [2/2]

template<class basis >
void AssembleDG< basis >::reinit ( const typename dealii::DoFHandler< lifex::dim >::active_cell_iterator new_cell,
const unsigned int  new_edge 
)

Reinitialize object on the current new_edge of the new_cell.

Definition at line 328 of file assemble_DG.hpp.

Member Data Documentation

◆ basis_ptr

template<class basis >
const std::unique_ptr<basis> AssembleDG< basis >::basis_ptr
protected

Basis function class.

Definition at line 72 of file assemble_DG.hpp.

◆ cell

template<class basis >
dealii::DoFHandler<lifex::dim>::active_cell_iterator AssembleDG< basis >::cell
protected

Cell object.

Definition at line 99 of file assemble_DG.hpp.

◆ dofs_per_cell

template<class basis >
unsigned int AssembleDG< basis >::dofs_per_cell
protected

Degree of freedom for each cell.

Definition at line 84 of file assemble_DG.hpp.

◆ face_handler

template<class basis >
FaceHandlerDG<lifex::dim> AssembleDG< basis >::face_handler
protected

Face element handler.

Definition at line 78 of file assemble_DG.hpp.

◆ face_handler_neigh

template<class basis >
FaceHandlerDG<lifex::dim> AssembleDG< basis >::face_handler_neigh
protected

Face element handler of the neighbor face.

Definition at line 81 of file assemble_DG.hpp.

◆ n_quad_points

template<class basis >
const unsigned int AssembleDG< basis >::n_quad_points
protected

Number of quadrature points in the volume: (n_quad_points_1D)^(dim).

Definition at line 90 of file assemble_DG.hpp.

◆ n_quad_points_1D

template<class basis >
const unsigned int AssembleDG< basis >::n_quad_points_1D
protected

Number of quadrature points in 1D.

Definition at line 87 of file assemble_DG.hpp.

◆ n_quad_points_face

template<class basis >
const unsigned int AssembleDG< basis >::n_quad_points_face
protected

Number of quadrature points in the face: (n_quad_points_1D)^(dim-1).

Definition at line 93 of file assemble_DG.hpp.

◆ poly_degree

template<class basis >
const unsigned int AssembleDG< basis >::poly_degree
protected

Polynomial degree.

Definition at line 96 of file assemble_DG.hpp.

◆ vol_handler

template<class basis >
VolumeHandlerDG<lifex::dim> AssembleDG< basis >::vol_handler
protected

Volume element handler.

Definition at line 75 of file assemble_DG.hpp.


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