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

Class to solve the monodomain equation with Fitzhugh-Nagumo ionic model for the cardiac electrophysiology using the discontinuous Galerkin method. More...

#include <monodomain_fhn_dg.hpp>

Inheritance diagram for DUBeat::models::MonodomainFHNDG< basis >:
Collaboration diagram for DUBeat::models::MonodomainFHNDG< basis >:

Public Member Functions

 MonodomainFHNDG ()
 Constructor. More...
 
- Public Member Functions inherited from ModelDG_t< basis >
 ModelDG_t (std::string model_name)
 Constructor. More...
 
 ModelDG_t (ModelDG_t< basis > &ModelDG_t)=default
 Default copy constructor. More...
 
 ModelDG_t (const ModelDG_t< basis > &ModelDG_t)=default
 Default const copy constructor. More...
 
 ModelDG_t (ModelDG_t< basis > &&ModelDG_t)=default
 Default move constructor. More...
 
virtual ~ModelDG_t ()=default
 Destructor. More...
 
- Public Member Functions inherited from ModelDG< basis >
 ModelDG (std::string model_name)
 Constructor. More...
 
 ModelDG (ModelDG< basis > &ModelDG)=default
 Default copy constructor. More...
 
 ModelDG (const ModelDG< basis > &ModelDG)=default
 Default const copy constructor. More...
 
 ModelDG (ModelDG< basis > &&ModelDG)=default
 Default move constructor. More...
 
unsigned int get_dofs_per_cell () const
 Return the number of degrees of freedom per element. More...
 
virtual ~ModelDG ()=default
 Destructor. More...
 

Private Member Functions

virtual void declare_parameters (lifex::ParamHandler &params) const override
 Override for declaration of additional parameters. More...
 
virtual void parse_parameters (lifex::ParamHandler &params) override
 Override to parse additional parameters. More...
 
void run () override
 Override for the simulation run. More...
 
void update_time () override
 Override to update time for both u and w. More...
 
void time_initialization () override
 Override to initialize both u and w. More...
 
void assemble_system () override
 Assembly of the Monodomain system. More...
 

Private Attributes

double ChiM
 Monodomain equation parameter. More...
 
double Sigma
 Diffusion scalar parameter. More...
 
double Cm
 Membrane capacity. More...
 
double kappa
 Factor for the nonlinear reaction in Fitzhugh Nagumo model. More...
 
double epsilon
 ODE parameter. More...
 
double gamma
 ODE parameter. More...
 
double a
 ODE parameter. More...
 
lifex::LinAlg::MPI::Vector solution_owned_w
 Solution gating variable, without ghost entries. More...
 
lifex::LinAlg::MPI::Vector solution_w
 Solution gating variable, with ghost entries. More...
 
lifex::LinAlg::MPI::Vector solution_ex_owned_w
 Solution exact gating variable, without ghost entries. More...
 
lifex::LinAlg::MPI::Vector solution_ex_w
 Solution exact gating variable, without ghost entries. More...
 
std::shared_ptr< lifex::utils::FunctionDirichlet > w_ex
 dealii::Pointer to exact solution function gating variable. More...
 
std::shared_ptr< dealii::Function< lifex::dim > > grad_w_ex
 dealii::Pointer to exact gradient solution Function gating variable More...
 
lifex::utils::BDFHandler< lifex::LinAlg::MPI::Vector > bdf_handler_w
 BDF time advancing handler. More...
 
lifex::LinAlg::MPI::Vector solution_bdf_w
 BDF solution, with ghost entries. More...
 
lifex::LinAlg::MPI::Vector solution_ext_w
 BDF extrapolated solution, with ghost entries. More...
 
lifex::LinAlg::MPI::SparseMatrix matrix_t0
 Component of the system matrix that does not depend on time. More...
 

Additional Inherited Members

- Protected Member Functions inherited from ModelDG_t< basis >
virtual void intermediate_error_print (const lifex::LinAlg::MPI::Vector &solution_owned, const lifex::LinAlg::MPI::Vector &solution_ex_owned, const std::shared_ptr< dealii::Function< lifex::dim >> &u_ex, const char *solution_name=(char *)"u")
 Computation of the \(L^\infty\) error at an intermediate time-step. More...
 
void run () override
 Override for the simulation run. More...
 
- Protected Member Functions inherited from ModelDG< basis >
virtual void setup_system ()
 Setup of the problem before the resolution. More...
 
void make_sparsity_pattern (const DoFHandlerDG< basis > &dof, dealii::DynamicSparsityPattern &sparsity, const dealii::AffineConstraints< double > &constraints=dealii::AffineConstraints< double >(), const bool keep_constrained_dofs=true, const dealii::types::subdomain_id subdomain_id=dealii::numbers::invalid_subdomain_id)
 Creation of the sparsity pattern to assign to the system matrix before assembling. More...
 
void initialize_solution (lifex::LinAlg::MPI::Vector &solution_owned, lifex::LinAlg::MPI::Vector &solution)
 To inizialize the solutions using the deal.II reinit. More...
 
void create_mesh ()
 Load the mesh from the default path. More...
 
void create_mesh (std::string mesh_path)
 Load the mesh from a user-defined path. More...
 
void solve_system ()
 System solving. More...
 
void compute_errors (const lifex::LinAlg::MPI::Vector &solution_owned, const lifex::LinAlg::MPI::Vector &solution_ex_owned, const std::shared_ptr< dealii::Function< lifex::dim >> &u_ex, const std::shared_ptr< dealii::Function< lifex::dim >> &grad_u_ex, const char *solution_name=(char *)"u") const
 To compute the \(L^\infty\) error, the \(L^2\) error, the \(H^1\) error and the \(DG\) error at the end of system solving, it exploits the DGComputeErrors<basis> class. More...
 
virtual void output_results (std::string output_name="solution") const
 Output of results. More...
 
virtual void conversion_to_fem (lifex::LinAlg::MPI::Vector &sol_owned)
 To convert a discretized solution from modal to nodal basis (does nothing if problem is already in nodal basis), in-place version. More...
 
virtual lifex::LinAlg::MPI::Vector conversion_to_fem (const lifex::LinAlg::MPI::Vector &sol_owned) const
 To convert a discretized solution from modal to nodal basis (does nothing if problem is already in nodal basis), const version. More...
 
virtual void conversion_to_fem (lifex::LinAlg::MPI::Vector &sol_owned, const std::string fem_file_path, const unsigned int degree_fem=1, const double scaling_factor=1) const
 
virtual void conversion_to_dub (lifex::LinAlg::MPI::Vector &sol_owned)
 To convert a discretized solution in Dubiner basis (only for problems using Dubiner basis). More...
 
virtual lifex::LinAlg::MPI::Vector conversion_to_dub (const lifex::LinAlg::MPI::Vector &sol_owned) const
 To convert a discretized solution in Dubiner basis (only for problems using Dubiner basis). More...
 
virtual void conversion_to_dub (lifex::LinAlg::MPI::Vector &sol_owned, const std::string fem_file_path, const unsigned int degree_fem=1, const double scaling_factor=1) const
 
void discretize_analytical_solution (const std::shared_ptr< dealii::Function< lifex::dim >> &u_analytical, lifex::LinAlg::MPI::Vector &sol_owned)
 Conversion of an analytical solution from FEM to basis coefficients. More...
 
void conversion_to_fem (lifex::LinAlg::MPI::Vector &sol_owned)
 Conversion of a discretized solution from Dubiner coefficients to FEM coefficients. More...
 
lifex::LinAlg::MPI::Vector conversion_to_fem (const lifex::LinAlg::MPI::Vector &sol_owned) const
 Conversion to FEM coefficients, const version. More...
 
void conversion_to_fem (lifex::LinAlg::MPI::Vector &sol_owned, const std::string fem_file_path, const unsigned int degree_fem, const double scaling_factor) const
 
void conversion_to_dub (lifex::LinAlg::MPI::Vector &sol_owned)
 Conversion of a discretized solution from FEM coefficients to Dubiner coefficients. More...
 
lifex::LinAlg::MPI::Vector conversion_to_dub (const lifex::LinAlg::MPI::Vector &sol_owned) const
 Conversion to DUB coefficients, const version. More...
 
void conversion_to_dub (lifex::LinAlg::MPI::Vector &sol_owned, const std::string fem_file_path, const unsigned int degree_fem, const double scaling_factor) const
 
void discretize_analytical_solution (const std::shared_ptr< dealii::Function< lifex::dim >> &u_analytical, lifex::LinAlg::MPI::Vector &sol_owned)
 Conversion of an analytical solution from FEM to basis coefficients. More...
 
void discretize_analytical_solution (const std::shared_ptr< dealii::Function< lifex::dim >> &u_analytical, lifex::LinAlg::MPI::Vector &sol_owned)
 Conversion of an analytical solution from FEM to basis coefficients. More...
 
- Protected Attributes inherited from ModelDG_t< basis >
double prm_time_init
 Initial time. More...
 
double prm_time_final
 Final time. More...
 
double prm_time_step
 Time-step amplitude. More...
 
unsigned int prm_bdf_order
 BDF order. More...
 
double time
 Current time. More...
 
unsigned int timestep_number
 Current time-step number. More...
 
lifex::utils::BDFHandler< lifex::LinAlg::MPI::Vector > bdf_handler
 BDF time advancing handler. More...
 
lifex::LinAlg::MPI::Vector solution_bdf
 BDF solution, with ghost entries. More...
 
lifex::LinAlg::MPI::Vector solution_ext
 BDF extrapolated solution, with ghost entries. More...
 
- Protected Attributes inherited from ModelDG< basis >
const std::string model_name
 Name of the class/problem. More...
 
unsigned int prm_fe_degree
 Polynomials degree. More...
 
unsigned int prm_n_refinements
 Mesh refinement level (>=1). More...
 
double prm_penalty_coeff
 DG Penalty coefficient. More...
 
double prm_stability_coeff
 DG stabilty coefficient. More...
 
double scaling_factor = 1
 Scaling factor. More...
 
std::shared_ptr< lifex::utils::MeshHandler > triangulation
 Triangulation (internal use for useful already implemented methods). More...
 
unsigned int dofs_per_cell
 Number of degrees of freedom per cell. More...
 
DoFHandlerDG< basis > dof_handler
 DoFHandler (internal use for useful already implemented methods). More...
 
std::shared_ptr< DUBFEMHandler< basis > > dub_fem_values
 Member used for conversions between analytical, nodal and modal representations of the solutions. More...
 
std::unique_ptr< AssembleDG< basis > > assemble
 Matrix assembler. More...
 
lifex::utils::LinearSolverHandler< lifex::LinAlg::MPI::Vector > linear_solver
 Linear solver handler. More...
 
lifex::utils::PreconditionerHandler preconditioner
 Preconditioner handler. More...
 
lifex::LinAlg::MPI::SparseMatrix matrix
 Distributed matrix of the linear system. More...
 
lifex::LinAlg::MPI::Vector rhs
 Distributed right hand side vector of the linear system. More...
 
lifex::LinAlg::MPI::Vector solution_owned
 Distributed solution vector, without ghost entries. More...
 
lifex::LinAlg::MPI::Vector solution
 Distributed solution vector, with ghost entries. More...
 
lifex::LinAlg::MPI::Vector solution_ex_owned
 Distributed exact solution vector, without ghost entries. More...
 
lifex::LinAlg::MPI::Vector solution_ex
 Distributed exact solution vector, without ghost entries. More...
 
std::shared_ptr< lifex::utils::FunctionDirichlet > u_ex
 Pointer to exact solution function. More...
 
std::shared_ptr< dealii::Function< lifex::dim > > grad_u_ex
 Pointer to exact gradient solution Function. More...
 
std::shared_ptr< dealii::Function< lifex::dim > > f_ex
 Known forcing term. More...
 
std::shared_ptr< dealii::Function< lifex::dim > > g_n
 Neumann boundary conditions. More...
 

Detailed Description

template<class basis>
class DUBeat::models::MonodomainFHNDG< basis >

Class to solve the monodomain equation with Fitzhugh-Nagumo ionic model for the cardiac electrophysiology using the discontinuous Galerkin method.

\[ \begin{aligned} \chi_m C_m \frac{\partial V_m}{\partial t} -\nabla \cdot (\Sigma \nabla V_m)+\chi_m I_{ion}(V_m, w) &= I_{ext} & \quad & \text{in } \Omega = (-1, 1)^d \times (0, T], \\ I_{ion}(V_m, w) &= k V_m (V_m - a)(V_m - 1) + w & \quad & \text{in } \Omega \times (0, T], \\ \frac{\partial w}{\partial t} &= \epsilon (V_m - \gamma w) & \quad & \text{in } \Omega \times (0, T], \\ \Sigma \frac{\partial V_m}{\partial n} &= g & \quad & \text{on } \partial\Omega \times (0, T], \\ V_m &= V_{m_\mathrm{ex}} & \quad & \text{in } \Omega \times \{ t = 0 \}, \\ w &= w_\mathrm{ex} & \quad & \text{in } \Omega \times \{ t = 0 \} \end{aligned} \]

where \( I_{ion}(V_m, w)\) is defined through the FitzHugh-Nagumo model.

In particular, it can be solved using the Lagrangian basis (basis=dealii::FE_SimplexDGP<lifex::dim>) or the Dubiner basis (basis=DUBValues<lifex::dim>).

The problem is time-discretized using the implicit finite difference scheme \(\mathrm{BDF}\sigma\) (where \(\sigma = 1,2,...\) is the order of the BDF formula) as follows:

\[ \begin{aligned} \chi_m C_m ( \frac{\alpha_{\mathrm{BDF}\sigma} V_m^{n+1} - V_{m_{\mathrm{BDF}\sigma}}^n}{\Delta t}) -\nabla \cdot (\Sigma \nabla V_m^{n+1})+\chi_m I_{ion}(V_m^{n+1}, w^n) &= I_{ext}^{n+1} & \quad & \text{in } \Omega = (-1, 1)^d \times \{n=0, \dots, N \}, \\ I_{ion}(V_m^{n+1}, w^n) &= k V_m^{n+1} (V_m^n - a)(V_m^n - 1) + w^n & \quad & \text{in } \Omega \times \{n=0, \dots, N \}, \\ \frac{\alpha_{\mathrm{BDF}\sigma} w^{n+1} - w_{\mathrm{BDF}\sigma}^n}{\Delta t} &= \epsilon (V_m^{n+1} - \gamma w^n) & \quad & \text{in } \Omega \times \{n=0, \dots, N \}, \\ \Sigma \frac{\partial V_m^{n+1}}{\partial n} &= g & \quad & \text{on } \partial\Omega \times \{n=0, \dots, N \}, \\ V_m^0 &= V_{m_\mathrm{ex}} & \quad & \text{in } \Omega \times \{ n = 0 \}, \\ w^0 &= w_\mathrm{ex} & \quad & \text{in } \Omega \times \{ n = 0 \} \end{aligned} \]

where \(\Delta t = t^{n+1}-t^{n}\) is the time step.

Boundary conditions, initial conditions and source terms are provided assuming that the exact solution is:

\[ \begin{aligned} d=2: \: &V_{m_\mathrm{ex}}(x,y) &= \sin(2\pi x)\sin(2\pi y)e^{-5t}, \hspace{6mm} &(x,y) &\in \Omega=(1,1)^2&, t \in [0,T], \\ &w_{\mathrm{ex}}(x,y) &= \frac{\epsilon}{\epsilon\cdot\gamma -5}\sin(2\pi x)\sin(2\pi y)e^{-5t}, \hspace{6mm} &(x,y) &\in \Omega=(1,1)^2&, t \in [0,T], \\ d=3: \: &V_{m_\mathrm{ex}}(x,y,z) &= \sin\left(2\pi x + \frac{\pi}{4}\right)\sin\left(2\pi y + \frac{\pi}{4}\right)\sin\left(2\pi z + \frac{\pi}{4}\right) e^{-5t}, \hspace{6mm} &(x,y,z) &\in \Omega=(1,1)^3&, t \in [0,T], \\ &w_{\mathrm{ex}}(x,y,z) &= \frac{\epsilon}{\epsilon\cdot\gamma -5} \sin(2\pi x)\sin(2\pi y)\sin(2\pi z) e^{-5t}, \hspace{6mm} &(x,y,z) &\in \Omega=(1,1)^3&, t \in [0,T]. \end{aligned} \]

Finally, \(d\) is specified in the lifex configuration and \(T\) as well as the monodomain scalar parameters in the .prm parameter file.

Definition at line 446 of file monodomain_fhn_dg.hpp.

Constructor & Destructor Documentation

◆ MonodomainFHNDG()

template<class basis >
DUBeat::models::MonodomainFHNDG< basis >::MonodomainFHNDG ( )
inline

Constructor.

Definition at line 340 of file monodomain_fhn_dg.hpp.

Member Function Documentation

◆ assemble_system()

template<class basis >
void DUBeat::models::MonodomainFHNDG< basis >::assemble_system
overrideprivatevirtual

Assembly of the Monodomain system.

Implements ModelDG< basis >.

Definition at line 780 of file monodomain_fhn_dg.hpp.

◆ declare_parameters()

template<class basis >
void DUBeat::models::MonodomainFHNDG< basis >::declare_parameters ( lifex::ParamHandler &  params) const
overrideprivatevirtual

Override for declaration of additional parameters.

Reimplemented from ModelDG_t< basis >.

Definition at line 518 of file monodomain_fhn_dg.hpp.

◆ parse_parameters()

template<class basis >
void DUBeat::models::MonodomainFHNDG< basis >::parse_parameters ( lifex::ParamHandler &  params)
overrideprivatevirtual

Override to parse additional parameters.

Reimplemented from ModelDG_t< basis >.

Definition at line 612 of file monodomain_fhn_dg.hpp.

◆ run()

template<class basis >
void DUBeat::models::MonodomainFHNDG< basis >::run
overrideprivatevirtual

Override for the simulation run.

Reimplemented from ModelDG< basis >.

Definition at line 662 of file monodomain_fhn_dg.hpp.

◆ time_initialization()

template<class basis >
void DUBeat::models::MonodomainFHNDG< basis >::time_initialization
overrideprivatevirtual

Override to initialize both u and w.

Reimplemented from ModelDG_t< basis >.

Definition at line 746 of file monodomain_fhn_dg.hpp.

◆ update_time()

template<class basis >
void DUBeat::models::MonodomainFHNDG< basis >::update_time
overrideprivatevirtual

Override to update time for both u and w.

Reimplemented from ModelDG_t< basis >.

Definition at line 724 of file monodomain_fhn_dg.hpp.

Member Data Documentation

◆ a

template<class basis >
double DUBeat::models::MonodomainFHNDG< basis >::a
private

ODE parameter.

Definition at line 468 of file monodomain_fhn_dg.hpp.

◆ bdf_handler_w

template<class basis >
lifex::utils::BDFHandler<lifex::LinAlg::MPI::Vector> DUBeat::models::MonodomainFHNDG< basis >::bdf_handler_w
private

BDF time advancing handler.

Definition at line 482 of file monodomain_fhn_dg.hpp.

◆ ChiM

template<class basis >
double DUBeat::models::MonodomainFHNDG< basis >::ChiM
private

Monodomain equation parameter.

Definition at line 456 of file monodomain_fhn_dg.hpp.

◆ Cm

template<class basis >
double DUBeat::models::MonodomainFHNDG< basis >::Cm
private

Membrane capacity.

Definition at line 460 of file monodomain_fhn_dg.hpp.

◆ epsilon

template<class basis >
double DUBeat::models::MonodomainFHNDG< basis >::epsilon
private

ODE parameter.

Definition at line 464 of file monodomain_fhn_dg.hpp.

◆ gamma

template<class basis >
double DUBeat::models::MonodomainFHNDG< basis >::gamma
private

ODE parameter.

Definition at line 466 of file monodomain_fhn_dg.hpp.

◆ grad_w_ex

template<class basis >
std::shared_ptr<dealii::Function<lifex::dim> > DUBeat::models::MonodomainFHNDG< basis >::grad_w_ex
private

dealii::Pointer to exact gradient solution Function gating variable

Definition at line 480 of file monodomain_fhn_dg.hpp.

◆ kappa

template<class basis >
double DUBeat::models::MonodomainFHNDG< basis >::kappa
private

Factor for the nonlinear reaction in Fitzhugh Nagumo model.

Definition at line 462 of file monodomain_fhn_dg.hpp.

◆ matrix_t0

template<class basis >
lifex::LinAlg::MPI::SparseMatrix DUBeat::models::MonodomainFHNDG< basis >::matrix_t0
private

Component of the system matrix that does not depend on time.

Definition at line 489 of file monodomain_fhn_dg.hpp.

◆ Sigma

template<class basis >
double DUBeat::models::MonodomainFHNDG< basis >::Sigma
private

Diffusion scalar parameter.

Definition at line 458 of file monodomain_fhn_dg.hpp.

◆ solution_bdf_w

template<class basis >
lifex::LinAlg::MPI::Vector DUBeat::models::MonodomainFHNDG< basis >::solution_bdf_w
private

BDF solution, with ghost entries.

Definition at line 484 of file monodomain_fhn_dg.hpp.

◆ solution_ex_owned_w

template<class basis >
lifex::LinAlg::MPI::Vector DUBeat::models::MonodomainFHNDG< basis >::solution_ex_owned_w
private

Solution exact gating variable, without ghost entries.

Definition at line 474 of file monodomain_fhn_dg.hpp.

◆ solution_ex_w

template<class basis >
lifex::LinAlg::MPI::Vector DUBeat::models::MonodomainFHNDG< basis >::solution_ex_w
private

Solution exact gating variable, without ghost entries.

Definition at line 476 of file monodomain_fhn_dg.hpp.

◆ solution_ext_w

template<class basis >
lifex::LinAlg::MPI::Vector DUBeat::models::MonodomainFHNDG< basis >::solution_ext_w
private

BDF extrapolated solution, with ghost entries.

Definition at line 486 of file monodomain_fhn_dg.hpp.

◆ solution_owned_w

template<class basis >
lifex::LinAlg::MPI::Vector DUBeat::models::MonodomainFHNDG< basis >::solution_owned_w
private

Solution gating variable, without ghost entries.

Definition at line 470 of file monodomain_fhn_dg.hpp.

◆ solution_w

template<class basis >
lifex::LinAlg::MPI::Vector DUBeat::models::MonodomainFHNDG< basis >::solution_w
private

Solution gating variable, with ghost entries.

Definition at line 472 of file monodomain_fhn_dg.hpp.

◆ w_ex

template<class basis >
std::shared_ptr<lifex::utils::FunctionDirichlet> DUBeat::models::MonodomainFHNDG< basis >::w_ex
private

dealii::Pointer to exact solution function gating variable.

Definition at line 478 of file monodomain_fhn_dg.hpp.


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