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

Class to compute the errors between the numerical solution (solution_owned) and the exact solution (solution_ex_owned). More...

#include <compute_errors_DG.hpp>

Public Member Functions

 ComputeErrorsDG (const unsigned int degree, const double stability_coeff, const unsigned int local_dofs, const DoFHandlerDG< basis > &dof_hand, const unsigned int nref_, const std::string &title)
 Constructor. More...
 
 ComputeErrorsDG (ComputeErrorsDG< basis > &ComputeErrorsDG)=default
 Default copy constructor. More...
 
 ComputeErrorsDG (const ComputeErrorsDG< basis > &ComputeErrorsDG)=default
 Default const copy constructor. More...
 
 ComputeErrorsDG (ComputeErrorsDG< basis > &&ComputeErrorsDG)=default
 Default move constructor. More...
 
void reinit (const lifex::LinAlg::MPI::Vector &sol_owned, const lifex::LinAlg::MPI::Vector &sol_ex_owned, const std::shared_ptr< dealii::Function< lifex::dim >> &u_ex_input, const std::shared_ptr< dealii::Function< lifex::dim >> &grad_u_ex_input, const char *solution_name_input)
 Reinitialization with new computed and exact solutions. More...
 
void compute_errors (std::list< const char * > errors_defs={"inf", "L2", "H1", "DG"})
 Compute errors following the preferences in the list. More...
 
std::vector< double > output_errors (std::list< const char * > errors_defs={"inf", "L2", "H1", "DG"}) const
 Output of the errors following the preferences in the list. More...
 
void update_datafile () const
 Update the datafile with the new errors. More...
 
void initialize_datafile () const
 Create the datafile for the errors. More...
 

Private Member Functions

void compute_error_inf ()
 Compute the \(L^\infty\) error. More...
 
void compute_error_L2 ()
 Compute the \(L^2\) error. More...
 
void compute_error_H1 ()
 Compute the \(H^1\) error. More...
 
void compute_error_DG ()
 Compute the \(DG\) error. More...
 
std::string get_date () const
 Return the current date and time. More...
 
void compute_error_inf ()
 Specialized version for Dubiner basis. More...
 

Private Attributes

const DoFHandlerDG< basis > & dof_handler
 Dof handler object of the problem. More...
 
const unsigned int n_quad_points
 Number of quadrature points in the volume element. More...
 
const unsigned int n_quad_points_face
 Number of quadrature points on the face element. More...
 
const unsigned int dofs_per_cell
 Number of degrees of freedom per cell. More...
 
const double stability_coefficient
 Stability coefficient, needed for the computation of the \(DG\) error. More...
 
const unsigned int fe_degree
 Polynomial degree. More...
 
const unsigned int nref
 Number of refinements. More...
 
const std::string model_name
 Name of the model. More...
 
const std::unique_ptr< basis > basis_ptr
 Pointer to the basis handler. More...
 
lifex::LinAlg::MPI::Vector solution_owned
 Computed solution. More...
 
lifex::LinAlg::MPI::Vector solution_ex_owned
 Exact solution to compare with the numerical one. More...
 
std::shared_ptr< dealii::Function< lifex::dim > > u_ex
 Analytical exact solution. More...
 
std::shared_ptr< dealii::Function< lifex::dim > > grad_u_ex
 Gradient of the analytical exact solution. More...
 
char * solution_name
 String that contains the solution name (to write on output files, default "u"). More...
 
std::array< double, 4 > errors
 Array that contains the current errors (in order, \(L^\infty\), \(L^2\), \(H^1\), \(DG\)). More...
 
std::shared_ptr< DUBFEMHandler< basis > > dub_fem_values
 Member to compute conversion to FEM basis, needed for L^inf error. More...
 

Detailed Description

template<class basis>
class ComputeErrorsDG< basis >

Class to compute the errors between the numerical solution (solution_owned) and the exact solution (solution_ex_owned).

Four definitions of error are implemented:

\[ \|u-u_h\|_{L^\infty(\Omega)}:= \sup_{x\in \Omega} |u(x)-u_h(x)| \]

\[ \|u-u_h\|_{L^2(\Omega)}^2 := \int_{\Omega} |u(x)-u_h(x)|^2 \, dx \]

\[ \|u-u_h\|_{H^1(\Omega)}^2 := \|u-u_h\|_{L^2(\Omega)}^2 + \|\nabla u-\nabla u_h\|_{L^2(\Omega)}^2 \]

\[ \|u-u_h\|_{DG(\Omega)}^2 := \|\nabla u-\nabla u_h\|_{L^2(\Omega)}^2 +\gamma \|[[ u- u_h]]\|_{L^2(\mathcal{F})}^2 \]

where \(\gamma \) is the stability coefficient and the \(L^2(\mathcal{F}) \) norm is computed on the faces instead of the volume. Furthermore, the class can create and update a datafile containing the errors from different simulations that has been previously run with different polynomial orders and grid refinements.

Definition at line 71 of file compute_errors_DG.hpp.

Constructor & Destructor Documentation

◆ ComputeErrorsDG() [1/4]

template<class basis >
ComputeErrorsDG< basis >::ComputeErrorsDG ( const unsigned int  degree,
const double  stability_coeff,
const unsigned int  local_dofs,
const DoFHandlerDG< basis > &  dof_hand,
const unsigned int  nref_,
const std::string &  title 
)
inline

Constructor.

Definition at line 569 of file compute_errors_DG.hpp.

◆ ComputeErrorsDG() [2/4]

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

Default copy constructor.

◆ ComputeErrorsDG() [3/4]

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

Default const copy constructor.

◆ ComputeErrorsDG() [4/4]

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

Default move constructor.

Member Function Documentation

◆ compute_error_DG()

template<class basis >
void ComputeErrorsDG< basis >::compute_error_DG
private

Compute the \(DG\) error.

Definition at line 441 of file compute_errors_DG.hpp.

◆ compute_error_H1()

template<class basis >
void ComputeErrorsDG< basis >::compute_error_H1
private

Compute the \(H^1\) error.

Definition at line 375 of file compute_errors_DG.hpp.

◆ compute_error_inf() [1/2]

template<class basis >
void ComputeErrorsDG< basis >::compute_error_inf
private

Compute the \(L^\infty\) error.

Definition at line 298 of file compute_errors_DG.hpp.

◆ compute_error_inf() [2/2]

void ComputeErrorsDG< DUBValues< lifex::dim > >::compute_error_inf ( )
private

Specialized version for Dubiner basis.

Before computing the \(L^\infty\) error, the vector solutions are transformed in terms of FEM coefficients.

Definition at line 310 of file compute_errors_DG.hpp.

◆ compute_error_L2()

template<class basis >
void ComputeErrorsDG< basis >::compute_error_L2
private

Compute the \(L^2\) error.

Definition at line 325 of file compute_errors_DG.hpp.

◆ compute_errors()

template<class basis >
void ComputeErrorsDG< basis >::compute_errors ( std::list< const char * >  errors_defs = {"inf", "L2", "H1", "DG"})

Compute errors following the preferences in the list.

E.g., errors_defs={"L2"} will only compute the \(L^2\) error.

Definition at line 227 of file compute_errors_DG.hpp.

◆ get_date()

template<class basis >
std::string ComputeErrorsDG< basis >::get_date
private

Return the current date and time.

Definition at line 556 of file compute_errors_DG.hpp.

◆ initialize_datafile()

template<class basis >
void ComputeErrorsDG< basis >::initialize_datafile

Create the datafile for the errors.

Definition at line 531 of file compute_errors_DG.hpp.

◆ output_errors()

template<class basis >
std::vector< double > ComputeErrorsDG< basis >::output_errors ( std::list< const char * >  errors_defs = {"inf", "L2", "H1", "DG"}) const

Output of the errors following the preferences in the list.

E.g., errors_defs={"L2"} will only output the \(L^2\) error. The output vector contains the errors in the order of the input list.

Definition at line 272 of file compute_errors_DG.hpp.

◆ reinit()

template<class basis >
void ComputeErrorsDG< basis >::reinit ( const lifex::LinAlg::MPI::Vector &  sol_owned,
const lifex::LinAlg::MPI::Vector &  sol_ex_owned,
const std::shared_ptr< dealii::Function< lifex::dim >> &  u_ex_input,
const std::shared_ptr< dealii::Function< lifex::dim >> &  grad_u_ex_input,
const char *  solution_name_input 
)

Reinitialization with new computed and exact solutions.

Definition at line 211 of file compute_errors_DG.hpp.

◆ update_datafile()

template<class basis >
void ComputeErrorsDG< basis >::update_datafile

Update the datafile with the new errors.

Definition at line 569 of file compute_errors_DG.hpp.

Member Data Documentation

◆ basis_ptr

template<class basis >
const std::unique_ptr<basis> ComputeErrorsDG< basis >::basis_ptr
private

Pointer to the basis handler.

Definition at line 183 of file compute_errors_DG.hpp.

◆ dof_handler

template<class basis >
const DoFHandlerDG<basis>& ComputeErrorsDG< basis >::dof_handler
private

Dof handler object of the problem.

Definition at line 157 of file compute_errors_DG.hpp.

◆ dofs_per_cell

template<class basis >
const unsigned int ComputeErrorsDG< basis >::dofs_per_cell
private

Number of degrees of freedom per cell.

Definition at line 168 of file compute_errors_DG.hpp.

◆ dub_fem_values

template<class basis >
std::shared_ptr<DUBFEMHandler<basis> > ComputeErrorsDG< basis >::dub_fem_values
private

Member to compute conversion to FEM basis, needed for L^inf error.

Definition at line 206 of file compute_errors_DG.hpp.

◆ errors

template<class basis >
std::array<double, 4> ComputeErrorsDG< basis >::errors
private

Array that contains the current errors (in order, \(L^\infty\), \(L^2\), \(H^1\), \(DG\)).

Definition at line 203 of file compute_errors_DG.hpp.

◆ fe_degree

template<class basis >
const unsigned int ComputeErrorsDG< basis >::fe_degree
private

Polynomial degree.

Definition at line 174 of file compute_errors_DG.hpp.

◆ grad_u_ex

template<class basis >
std::shared_ptr<dealii::Function<lifex::dim> > ComputeErrorsDG< basis >::grad_u_ex
private

Gradient of the analytical exact solution.

Definition at line 195 of file compute_errors_DG.hpp.

◆ model_name

template<class basis >
const std::string ComputeErrorsDG< basis >::model_name
private

Name of the model.

Definition at line 180 of file compute_errors_DG.hpp.

◆ n_quad_points

template<class basis >
const unsigned int ComputeErrorsDG< basis >::n_quad_points
private

Number of quadrature points in the volume element.

By default: \((degree+2)^{dim}\).

Definition at line 161 of file compute_errors_DG.hpp.

◆ n_quad_points_face

template<class basis >
const unsigned int ComputeErrorsDG< basis >::n_quad_points_face
private

Number of quadrature points on the face element.

By default: \((degree+2)^{dim-1}\).

Definition at line 165 of file compute_errors_DG.hpp.

◆ nref

template<class basis >
const unsigned int ComputeErrorsDG< basis >::nref
private

Number of refinements.

Definition at line 177 of file compute_errors_DG.hpp.

◆ solution_ex_owned

template<class basis >
lifex::LinAlg::MPI::Vector ComputeErrorsDG< basis >::solution_ex_owned
private

Exact solution to compare with the numerical one.

Definition at line 189 of file compute_errors_DG.hpp.

◆ solution_name

template<class basis >
char* ComputeErrorsDG< basis >::solution_name
private

String that contains the solution name (to write on output files, default "u").

Definition at line 199 of file compute_errors_DG.hpp.

◆ solution_owned

template<class basis >
lifex::LinAlg::MPI::Vector ComputeErrorsDG< basis >::solution_owned
private

Computed solution.

Definition at line 186 of file compute_errors_DG.hpp.

◆ stability_coefficient

template<class basis >
const double ComputeErrorsDG< basis >::stability_coefficient
private

Stability coefficient, needed for the computation of the \(DG\) error.

Definition at line 171 of file compute_errors_DG.hpp.

◆ u_ex

template<class basis >
std::shared_ptr<dealii::Function<lifex::dim> > ComputeErrorsDG< basis >::u_ex
private

Analytical exact solution.

Definition at line 192 of file compute_errors_DG.hpp.


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