LifeV
IonicNoblePurkinje.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 IonicNoblePurkinje
29  @brief Ionic model of Noble for Purkinje cells
30 
31  Model as in
32  Keener, James, and James Sneyd.
33  Mathematical Physiology: I: Cellular Physiology. Vol. 1.
34  Springer, 2010.
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 _IONICNOBLEPURKINJE_H_
46 #define _IONICNOBLEPURKINJE_H_
47 
48 #include <lifev/electrophysiology/solver/IonicModels/ElectroIonicModel.hpp>
49 
50 #include <Teuchos_RCP.hpp>
51 #include <Teuchos_ParameterList.hpp>
52 #include "Teuchos_XMLParameterListHelpers.hpp"
53 
54 namespace LifeV
55 {
56 //! IonicModel - This class implements an ionic model.
57 
58 class IonicNoblePurkinje : 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 std::shared_ptr<VectorElemental> elvecPtr_Type;
67  typedef RegionMesh<LinearTetra> mesh_Type;
68  //@}
69 
70 
71 
72  //! @name Constructors & Destructor
73  //@{
74 
75  //! Constructor
76  IonicNoblePurkinje();
77 
78  /*!
79  * @param Epetra communicator
80  * @param list of parameters in an xml file
81  */
82  IonicNoblePurkinje ( Teuchos::ParameterList& parameterList );
83 
84  /*!
85  * @param IonicNoblePurkinje object
86  */
87  IonicNoblePurkinje ( const IonicNoblePurkinje& model );
88  //! Destructor
89  virtual ~IonicNoblePurkinje() {}
90 
91  //@}
92 
93  //! @name Overloads
94  //@{
95 
96  IonicNoblePurkinje& operator= ( const IonicNoblePurkinje& model );
97 
98  //@}
99 
100  //! @name Setters and getters
101  //@{
102 
103  //parameters getters and setters
104  inline const Real& gi() const
105  {
106  return M_gi;
107  }
108 
109 
110  inline const Real& Cm() const
111  {
112  return M_Cm;
113  }
114 
115  inline const Real& vK() const
116  {
117  return M_vK;
118  }
119 
120  inline const Real& vNa() const
121  {
122  return M_vNa;
123  }
124 
125  inline const Real& Itotal() const
126  {
127  return M_Itot;
128  }
129 
130 
131  inline void setgi ( const Real& p )
132  {
133  this->M_gi = p;
134  }
135 
136 
137  inline void setCm ( const Real& p )
138  {
139  this->M_Cm = p;
140  }
141 
142  inline void setvNa ( const Real& p )
143  {
144  this->M_vNa = p;
145  }
146 
147  inline void setvK ( const Real& p )
148  {
149  this->M_vK = p;
150  }
151 
152 
153  inline static Real GeneralFunctionAlphaAndBeta (const Real& V, const Real& C1, const Real& C2, const Real& C3, const Real& C4, const Real& C5, const Real& V_0)
154  {
155  if (! (C1 != 0 && C2 == 0) )
156  {
157  return (C1 * std::exp ( (V - V_0) / C2) + C3 * (V - V_0) ) / (1 + C4 * std::exp ( (V - V_0) / C5) ) ;
158  }
159  else
160  {
161  return (C1 + C3 * (V - V_0) ) / (1 + C4 * std::exp ( (V - V_0) / C5) ) ;
162  }
163  }
164 
165 
166  inline static Real mInf (const Real& V)
167  {
168  Real alpham = GeneralFunctionAlphaAndBeta (V, 0, 1, 0.1, -1, -15, -48);
169  Real betam = GeneralFunctionAlphaAndBeta (V, 0, 1, -0.12, -1, 5, -8);
170  Real taum = alpham + betam;
171  return alpham / (taum);
172  }
173 
174  inline static Real hInf (const Real& V)
175  {
176  Real alphah = GeneralFunctionAlphaAndBeta (V, 0.17, -20, 0, 0, 1, -90);
177  Real betah = GeneralFunctionAlphaAndBeta (V, 1, 0, 0, 1, -10, -42);
178  Real tauh = alphah + betah;
179  return alphah / (tauh);
180  }
181 
182  inline static Real nInf (const Real& V)
183  {
184  Real alphan = GeneralFunctionAlphaAndBeta (V, 0, 1, 0.0001, -1, -10, -50);
185  Real betan = GeneralFunctionAlphaAndBeta (V, 0.002, -80, 0, 0, 1, -90);
186  Real taun = alphan + betan;
187  return alphan / (taun);
188  }
189  //inline const short int& Size() const { return M_numberOfEquations; }
190  //@}
191 
192 
193  //Compute the rhs on a single node or for the 0D case
194  void computeGatingRhs ( const std::vector<Real>& v, std::vector<Real>& rhs);
195 
196  void computeRhs ( const std::vector<Real>& v, std::vector<Real>& rhs);
197 
198  // compute the rhs with state variable interpolation
199  Real computeLocalPotentialRhs ( const std::vector<Real>& v );
200 
201  void computeGatingVariablesWithRushLarsen ( std::vector<Real>& v, const Real dt );
202 
203 
204  //! Display information about the model
205  void showMe();
206 
207  //! Solves the ionic model
208  //virtual void solveXbModel( const vector_Type& Calcium,
209  // const Real timeStep )=0;
210  //@}
211 
212 private:
213  //! Model Parameters
214 
215  //! Chemical kinetics parameters
216  Real M_gi;
217  Real M_vNa;
218  Real M_vK;
219  Real M_Cm;
220  Real M_Itot;
221 
222 
223  //! Xb states == equivalent to the number of equations
224  //short int M_numberOfEquations;
225 
226  //@}
227 
228 }; // class IonicMinimalModel
229 
230 
231 
232 }
233 
234 #endif
void updateInverseJacobian(const UInt &iQuadPt)
double Real
Generic real data.
Definition: LifeV.hpp:175