LifeV
IonicFitzHughNagumo.hpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
7 
8  This file is part of LifeV.
9 
10  LifeV is free software; you can redistribute it and/or modify
11  it under the terms of the GNU Lesser General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  LifeV is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  Lesser General Public License for more details.
19 
20  You should have received a copy of the GNU Lesser General Public License
21  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 
27 /*!
28  @file IonicFitzHughNagumo
29  @brief Ionic model of FitzHugh-Nagumo
30 
31  Model taken as in
32  Franzone, Piero Colli, et al.
33  "Adaptivity in space and time for reaction-diffusion systems in electrocardiology."
34  SIAM Journal on Scientific Computing 28.3 (2006): 942-962.
35 
36 
37  @date 01-2013
38  @author Giacomo Rosilho de Souza <giacomo.rosilhodesouza@epfl.ch>
39 
40  @contributors
41  @mantainer Giacomo Rosilho de Souza <giacomo.rosilhodesouza@epfl.ch>
42  @last update 01-2013
43  */
44 
45 
46 #ifndef _IONICFITZHUGHNAGUMO_H_
47 #define _IONICFITZHUGHNAGUMO_H_
48 
49 #include <lifev/electrophysiology/solver/IonicModels/ElectroIonicModel.hpp>
50 
51 #include <lifev/core/array/MatrixSmall.hpp>
52 #include <lifev/core/array/VectorSmall.hpp>
53 #include <lifev/core/array/MatrixEpetra.hpp>
54 #include <lifev/core/array/VectorEpetra.hpp>
55 
56 #include <Teuchos_RCP.hpp>
57 #include <Teuchos_ParameterList.hpp>
58 #include "Teuchos_XMLParameterListHelpers.hpp"
59 
60 namespace LifeV
61 {
62 //! IonicModel - This class implements an ionic model.
63 
64 //! @name Type definitions
65 //@{
66 typedef ElectroIonicModel super;
71 //@}
72 
73 class IonicFitzHughNagumo : public virtual ElectroIonicModel
74 {
75 
76 public:
77  //! @name Constructors & Destructor
78  //@{
79 
80  //! Constructor
82 
83  /*!
84  * @param Epetra communicator
85  * @param list of parameters in an xml file
86  */
87  IonicFitzHughNagumo ( Teuchos::ParameterList& parameterList );
88 
89  /*!
90  * @param IonicFitzHughNagumo object
91  */
93  //! Destructor
94  virtual ~IonicFitzHughNagumo() {}
95 
96  //@}
97 
98  //! @name Overloads
99  //@{
100 
102 
103  //@}
104 
105  //! @name Setters and getters
106  //@{
107 
108  //parameters getters and setters
109  inline const Real& G() const
110  {
111  return M_G;
112  }
113  inline const Real& Vth() const
114  {
115  return M_Vth;
116  }
117  inline const Real& Vp() const
118  {
119  return M_Vp;
120  }
121  inline const Real& Eta1() const
122  {
123  return M_Eta1;
124  }
125  inline const Real& Eta2() const
126  {
127  return M_Eta2;
128  }
129  inline const Real& Eta3() const
130  {
131  return M_Eta3;
132  }
133  inline const Real& Eta() const
134  {
135  return M_Eta;
136  }
137  inline const Real& Gamma() const
138  {
139  return M_Gamma;
140  }
141 
142  inline void setG ( const Real& G )
143  {
144  this->M_G = G;
145  }
146  inline void setVth ( const Real& Vth )
147  {
148  this->M_Vth = Vth;
149  }
150  inline void setVp ( const Real& Vp )
151  {
152  this->M_Vp = Vp;
153  this->M_Eta = M_Eta2 / M_Vp;
154  }
155  inline void setEta1 ( const Real& Eta1 )
156  {
157  this->M_Eta1 = Eta1;
158  }
159  inline void setEta2 ( const Real& Eta2 )
160  {
161  this->M_Eta2 = Eta2;
162  this->M_Eta = M_Eta2 / M_Vp;
163  this->M_Gamma = M_Eta2 * M_Eta3;
164  }
165  inline void setEta3 ( const Real& Eta3 )
166  {
167  this->M_Eta3 = Eta3;
168  this->M_Gamma = M_Eta2 * M_Eta3;
169  }
170 
171  void setup ( Teuchos::ParameterList& parameterList );
172  //Compute the rhs on a single node or for the 0D case
173  //this is to make visible the computeRhs defined in the base class, otherwise it is hidden by the ones defined here
174  using ElectroIonicModel::computeRhs;
175  void computeGatingRhs ( const std::vector<Real>& v, std::vector<Real>& rhs);
176  void computeRhs ( const std::vector<Real>& v, std::vector<Real>& rhs);
177  void computeRhs ( const VectorSmall<2>& v, VectorSmall<2>& rhs);
178 
179  // compute the rhs with state variable interpolation
180  Real computeLocalPotentialRhs ( const std::vector<Real>& v );
181 
182  //compute the Jacobian
183  matrix_Type getJac (const vector_Type& v, Real h = 1.0e-8);
184  std::vector< std::vector<Real> > getJac (const std::vector<Real>& v, Real h = 1.0e-8);
185  MatrixSmall<2, 2> getJac (const VectorSmall<2>& v, Real h = 1.0e-8);
186 
187  //! Display information about the model
188  void showMe();
189 
190 
191 private:
192  //! Model Parameters
193 
194  //! Chemical kinetics parameters
203 
204 
205  //@}
206 
207 }; // class IonicFitzHughNagumo
208 
209 
210 inline ElectroIonicModel* createIonicFitzHughNagumo()
211 {
212  return new IonicFitzHughNagumo();
213 }
214 
215 namespace
216 {
218 }
219 
220 }
221 
222 #endif
VectorEpetra - The Epetra Vector format Wrapper.
virtual ~IonicFitzHughNagumo()
Destructor.
void setEta1(const Real &Eta1)
IonicFitzHughNagumo(const IonicFitzHughNagumo &model)
void importFromHDF5(std::string const &fileName, std::string const &matrixName="matrix")
Read a matrix from a HDF5 (.h5) file.
ElectroIonicModel * createIonicFitzHughNagumo()
void updateInverseJacobian(const UInt &iQuadPt)
std::shared_ptr< vector_Type > vectorPtr_Type
void computeRhs(const VectorSmall< 2 > &v, VectorSmall< 2 > &rhs)
void showMe()
Display information about the model.
void setEta3(const Real &Eta3)
RegionMesh< LinearTetra > mesh_Type
void computeGatingRhs(const std::vector< Real > &v, std::vector< Real > &rhs)
A new ionic model must have the following methods ///.
void setEta2(const Real &Eta2)
IonicFitzHughNagumo & operator=(const IonicFitzHughNagumo &model)
Operator.
double Real
Generic real data.
Definition: LifeV.hpp:175
Real computeLocalPotentialRhs(const std::vector< Real > &v)
This methods contains the actual evaluation of the rhs of the voltage equation only (0D version) ...
MatrixEpetra< Real > matrix_Type
VectorEpetra vector_Type
IonicFitzHughNagumo(Teuchos::ParameterList &parameterList)
ElectroIonicModel super
MatrixSmall< 2, 2 > getJac(const VectorSmall< 2 > &v, Real h=1.0e-8)
void setup(Teuchos::ParameterList &parameterList)
Methods.
void computeRhs(const std::vector< Real > &v, std::vector< Real > &rhs)
This methods contains the actual evaluation of the rhs of all state variablesin the model (0D version...
std::vector< std::vector< Real > > getJac(const std::vector< Real > &v, Real h=1.0e-8)
This methods computes the Jacobian numerically.
matrix_Type getJac(const vector_Type &v, Real h=1.0e-8)
This methods computes the Jacobian numerically.