LifeV
IonicAlievPanfilov.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 IonicAlievPanfilov
29  @brief Ionic model of Aliev-Panfilov
30 
31  Model from:
32  Nash, Martyn P., and Alexander V. Panfilov.
33  "Electromechanical model of excitable tissue to study reentrant cardiac arrhythmias."
34  Progress in biophysics and molecular biology 85.2 (2004): 501-522.
35 
36  @date 01-2013
37  @author Simone Rossi <simone.rossi@epfl.ch>
38 
39  @contributors
40  @mantainer Simone Rossi <simone.rossi@epfl.ch>
41  @last update 01-2013
42  */
43 
44 
45 #ifndef _IONICALIEVPANFILOV_H_
46 #define _IONICALIEVPANFILOV_H_
47 
48 #include <lifev/electrophysiology/solver/IonicModels/ElectroIonicModel.hpp>
49 
50 
51 
52 
53 namespace LifeV
54 {
55 //! IonicModel - This class implements an ionic model.
56 
57 
58 class IonicAlievPanfilov : public virtual ElectroIonicModel
59 {
60 
61 public:
62  //! @name Type definitions
63  //@{
64  typedef ElectroIonicModel super;
65  typedef std::shared_ptr<VectorEpetra> vectorPtr_Type;
66  typedef RegionMesh<LinearTetra> mesh_Type;
67  //@}
68 
69 
70 
71  //! @name Constructors & Destructor
72  //@{
73 
74  //! Constructor
75  IonicAlievPanfilov();
76 
77  /*!
78  * @param Epetra communicator
79  * @param list of parameters in an xml file
80  */
81  IonicAlievPanfilov ( Teuchos::ParameterList& parameterList );
82 
83  /*!
84  * @param IonicAlievPanfilov object
85  */
86  IonicAlievPanfilov ( const IonicAlievPanfilov& model );
87  //! Destructor
88  virtual ~IonicAlievPanfilov() {}
89 
90  //@}
91 
92  //! @name Overloads
93  //@{
94 
95  IonicAlievPanfilov& operator= ( const IonicAlievPanfilov& model );
96 
97  //@}
98 
99  //! @name Setters and getters
100  //@{
101 
102  //parameters getters and setters
103  inline const Real& Mu1() const
104  {
105  return M_mu1;
106  }
107  inline const Real& Mu2() const
108  {
109  return M_mu2;
110  }
111  inline const Real& K() const
112  {
113  return M_k;
114  }
115  inline const Real& A() const
116  {
117  return M_a;
118  }
119  inline const Real& Epsilon() const
120  {
121  return M_epsilon;
122  }
123 
124  inline void setMu1 ( const Real& mu1 )
125  {
126  this->M_mu1 = mu1;
127  }
128  inline void setMu2 ( const Real& mu2 )
129  {
130  this->M_mu2 = mu2;
131  }
132  inline void setK ( const Real& k )
133  {
134  this->M_k = k;
135  }
136  inline void setA ( const Real& a )
137  {
138  this->M_a = a;
139  }
140  inline void setEpsilon ( const Real& epsilon )
141  {
142  this->M_epsilon = epsilon;
143  }
144 
145 
146  void setup ( Teuchos::ParameterList& parameterList );
147  //@}
148 
149  //! @name Methods
150  //@{
151 
152  //Compute the rhs on a single node or for the 0D case
153  void computeGatingRhs ( const std::vector<Real>& v, std::vector<Real>& rhs);
154 
155  void computeRhs ( const std::vector<Real>& v, std::vector<Real>& rhs);
156 
157  //Compute the rhs on a mesh/ 3D case
158  // void computeRhs( const std::vector<vectorPtr_Type>& v, std::vector<vectorPtr_Type>& rhs );
159  //
160  // void computeRhs( const std::vector<vectorPtr_Type>& v, const VectorEpetra& Iapp, std::vector<vectorPtr_Type>& rhs );
161  //
162 
163  // compute the rhs with state variable interpolation
164  Real computeLocalPotentialRhs ( const std::vector<Real>& v);
165 
166  // void computePotentialRhs( const std::vector<vectorPtr_Type>& v,
167  // const VectorEpetra& Iapp,
168  // std::vector<vectorPtr_Type>& rhs,
169  // FESpace<mesh_Type, MapEpetra>& uFESpace );
170 
171  //! Display information about the model
172  void showMe();
173 
174  //! Solves the ionic model
175  //virtual void solveXbModel( const vector_Type& Calcium,
176  // const Real timeStep )=0;
177  //@}
178 
179 private:
180  //! Model Parameters
181 
182  //! Chemical kinetics parameters
183  Real M_mu1;
184  Real M_mu2;
185  Real M_epsilon;
186  Real M_k;
187  Real M_a;
188 
189 
190  //@}
191 
192 }; // class IonicAlievPanfilov
193 
194 
195 
196 inline ElectroIonicModel* createIonicAlievPanfilov()
197 {
198  return new IonicAlievPanfilov();
199 }
200 
201 namespace
202 {
204 }
205 
206 }
207 
208 #endif
void updateInverseJacobian(const UInt &iQuadPt)
double Real
Generic real data.
Definition: LifeV.hpp:175