LifeV
IonicMinimalModel.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 IonicMinimalModel
29  @brief Ionic model of Bueno-Orovio et al.
30 
31  Model as in
32  Bueno-Orovio, Alfonso, Elizabeth M. Cherry, and Flavio H. Fenton.
33  "Minimal model for human ventricular action potentials in tissue."
34  Journal of theoretical biology 253.3 (2008): 544-560.
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 _IONICMINIMALMODEL_H_
46 #define _IONICMINIMALMODEL_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 IonicMinimalModel : 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  IonicMinimalModel();
77 
78  /*!
79  * @param Epetra communicator
80  * @param list of parameters in an xml file
81  */
82  IonicMinimalModel ( Teuchos::ParameterList& parameterList );
83 
84  /*!
85  * @param IonicMinimalModel object
86  */
87  IonicMinimalModel ( const IonicMinimalModel& model );
88  //! Destructor
89  virtual ~IonicMinimalModel() {}
90 
91  //@}
92 
93  //! @name Overloads
94  //@{
95 
96  IonicMinimalModel& operator= ( const IonicMinimalModel& model );
97 
98  //@}
99 
100  //! @name Setters and getters
101  //@{
102 
103  //parameters getters and setters
104  inline const Real& Uo() const
105  {
106  return M_uo;
107  }
108  inline const Real& Uu() const
109  {
110  return M_uu;
111  }
112  inline const Real& Tetav() const
113  {
114  return M_tetav;
115  }
116  inline const Real& Tetaw() const
117  {
118  return M_tetaw;
119  }
120  inline const Real& Tetavm() const
121  {
122  return M_tetavm;
123  }
124  inline const Real& Tetao() const
125  {
126  return M_tetao;
127  }
128  inline const Real& Tauv1() const
129  {
130  return M_tauv1;
131  }
132  inline const Real& Tauv2() const
133  {
134  return M_tauv2;
135  }
136  inline const Real& Tauvp() const
137  {
138  return M_tauvp;
139  }
140  inline const Real& Tauw1() const
141  {
142  return M_tauw1;
143  }
144  inline const Real& Tauw2() const
145  {
146  return M_tauw2;
147  }
148  inline const Real& Kw() const
149  {
150  return M_kw;
151  }
152  inline const Real& Uw() const
153  {
154  return M_uw;
155  }
156  inline const Real& Tauwp() const
157  {
158  return M_tauwp;
159  }
160  inline const Real& Taufi() const
161  {
162  return M_taufi;
163  }
164  inline const Real& Tauo1() const
165  {
166  return M_tauo1;
167  }
168  inline const Real& Tauo2() const
169  {
170  return M_tauo2;
171  }
172  inline const Real& Tauso1() const
173  {
174  return M_tauso1;
175  }
176  inline const Real& Tauso2() const
177  {
178  return M_tauso2;
179  }
180  inline const Real& Kso() const
181  {
182  return M_kso;
183  }
184  inline const Real& Uso() const
185  {
186  return M_uso;
187  }
188  inline const Real& Taus1() const
189  {
190  return M_taus1;
191  }
192  inline const Real& Taus2() const
193  {
194  return M_taus2;
195  }
196  inline const Real& Ks() const
197  {
198  return M_ks;
199  }
200  inline const Real& Us() const
201  {
202  return M_us;
203  }
204  inline const Real& Tausi() const
205  {
206  return M_tausi;
207  }
208  inline const Real& Tauinf() const
209  {
210  return M_tauwinf;
211  }
212  inline const Real& Winfstart() const
213  {
214  return M_winfstar;
215  }
216 
217 
218  inline void setUo ( const Real& p )
219  {
220  this->M_uo = p;
221  }
222  inline void setUu ( const Real& p )
223  {
224  this->M_uu = p;
225  }
226  inline void setTetav ( const Real& p )
227  {
228  this->M_tetav = p;
229  }
230  inline void setTetaw ( const Real& p )
231  {
232  this->M_tetaw = p;
233  }
234  inline void setTetavm ( const Real& p )
235  {
236  this->M_tetavm = p;
237  }
238  inline void setTetao ( const Real& p )
239  {
240  this->M_tetao = p;
241  }
242  inline void setTauv1 ( const Real& p )
243  {
244  this->M_tauv1 = p;
245  }
246  inline void setTauv2 ( const Real& p )
247  {
248  this->M_tauv2 = p;
249  }
250  inline void setTauvp ( const Real& p )
251  {
252  this->M_tauvp = p;
253  }
254  inline void setTauw1 ( const Real& p )
255  {
256  this->M_tauw1 = p;
257  }
258  inline void setTauw2 ( const Real& p )
259  {
260  this->M_tauw2 = p;
261  }
262  inline void setKw ( const Real& p )
263  {
264  this->M_kw = p;
265  }
266  inline void setUw ( const Real& p )
267  {
268  this->M_uw = p;
269  }
270  inline void setTauwp ( const Real& p )
271  {
272  this->M_tauwp = p;
273  }
274  inline void setTaufi ( const Real& p )
275  {
276  this->M_taufi = p;
277  }
278  inline void setTauo1 ( const Real& p )
279  {
280  this->M_tauo1 = p;
281  }
282  inline void setTauo2 ( const Real& p )
283  {
284  this->M_tauo2 = p;
285  }
286  inline void setTauso1 ( const Real& p )
287  {
288  this->M_tauso1 = p;
289  }
290  inline void setTauso2 ( const Real& p )
291  {
292  this->M_tauso2 = p;
293  }
294  inline void setKso ( const Real& p )
295  {
296  this->M_kso = p;
297  }
298  inline void setUso ( const Real& p )
299  {
300  this->M_uso = p;
301  }
302  inline void setTaus1 ( const Real& p )
303  {
304  this->M_taus1 = p;
305  }
306  inline void setTaus2 ( const Real& p )
307  {
308  this->M_taus2 = p;
309  }
310  inline void setKs ( const Real& p )
311  {
312  this->M_ks = p;
313  }
314  inline void setUs ( const Real& p )
315  {
316  this->M_us = p;
317  }
318  inline void setTausi ( const Real& p )
319  {
320  this->M_tausi = p;
321  }
322  inline void setTauinf ( const Real& p )
323  {
324  this->M_tauwinf = p;
325  }
326  inline void setWinfstart ( const Real& p )
327  {
328  this->M_winfstar = p;
329  }
330 
331 
332  //compute SVI with ETA (Yuppi Doo!!!)
333  void computePotentialRhsSVI ( const std::vector<vectorPtr_Type>& v,
334  std::vector<vectorPtr_Type>& rhs,
335  FESpace<mesh_Type, MapEpetra>& uFESpace);
336  //inline const short int& Size() const { return M_numberOfEquations; }
337  //@}
338 
339  //! @name Methods
340  //@{
341 
342  inline static Real Heaviside ( const Real& x)
343  {
344  if ( x > 0 )
345  {
346  return 1.0;
347  }
348  else
349  {
350  return 0.0;
351  }
352  }
353  //Compute the rhs on a single node or for the 0D case
354  void computeGatingRhs ( const std::vector<Real>& v, std::vector<Real>& rhs);
355 
356  void computeRhs ( const std::vector<Real>& v, std::vector<Real>& rhs);
357 
358  // compute the rhs with state variable interpolation
359  Real computeLocalPotentialRhs ( const std::vector<Real>& v );
360 
361  //
362  void computeGatingVariablesWithRushLarsen ( std::vector<Real>& v, const Real dt );
363 
364  //! Display information about the model
365  void showMe();
366 
367  //! Solves the ionic model
368  //virtual void solveXbModel( const vector_Type& Calcium,
369  // const Real timeStep )=0;
370  //@}
371 
372 private:
373  //! Model Parameters
374 
375  //! Chemical kinetics parameters
376  Real M_uo;
377  Real M_uu;
378  Real M_tetav;
379  Real M_tetaw;
380  Real M_tetavm;
381  Real M_tetao;
382  Real M_tauv1;
383  Real M_tauv2;
384  Real M_tauvp;
385  Real M_tauw1;
386  Real M_tauw2;
387  Real M_kw;
388  Real M_uw;
389  Real M_tauwp;
390  Real M_taufi;
391  Real M_tauo1;
392  Real M_tauo2;
393  Real M_tauso1;
394  Real M_tauso2;
395  Real M_kso;
396  Real M_uso;
397  Real M_taus1;
398  Real M_taus2;
399  Real M_ks;
400  Real M_us;
401  Real M_tausi;
402  Real M_tauwinf;
403  Real M_winfstar;
404 
405 
406  //! Xb states == equivalent to the number of equations
407  //short int M_numberOfEquations;
408 
409  //@}
410 
411 }; // class IonicMinimalModel
412 
413 
414 
415 inline ElectroIonicModel* createIonicMinimalModel()
416 {
417  return new IonicMinimalModel();
418 }
419 
420 namespace
421 {
423 }
424 
425 
426 
427 class MMTanhFunctor
428 {
429 public:
430  typedef Real return_Type;
431 
432  return_Type operator() (const VectorSmall<1>& p)
433  {
434  return std::tanh ( p[0] );
435  }
436  return_Type operator() (const Real& p)
437  {
438  return std::tanh ( p );
439  }
440 
441 
442  MMTanhFunctor() {}
443  MMTanhFunctor (const MMTanhFunctor&) {}
444  ~MMTanhFunctor() {}
445 };
446 
447 class MMHFunctor
448 {
449 public:
450  typedef Real return_Type;
451 
452  return_Type operator() (const VectorSmall<1>& p)
453  {
454  if (p[0] > 0.0)
455  {
456  return 1.0;
457  }
458  else
459  {
460  return 0.0;
461  }
462  }
463  return_Type operator() (const Real& p)
464  {
465  if (p > 0.0)
466  {
467  return 1.0;
468  }
469  else
470  {
471  return 0.0;
472  }
473  }
474 
475 
476  MMHFunctor() {}
477  MMHFunctor (const MMHFunctor&) {}
478  ~MMHFunctor() {}
479 };
480 
481 class MMSV
482 {
483 public:
484  typedef Real return_Type;
485 
486  return_Type operator() (const VectorSmall<1>& p)
487  {
488  std::cout.precision (15);
489  std::cout << "\nvalue: " << p[0];
490  }
491  return_Type operator() (const Real& p)
492  {
493  std::cout.precision (15);
494  std::cout << "\nvalue: " << p;
495  }
496 
497 
498  MMSV() {}
499  MMSV (const MMSV&) {}
500  ~MMSV() {}
501 };
502 
503 
504 
505 }
506 
507 #endif
void updateInverseJacobian(const UInt &iQuadPt)
Epetra_Import const & importer()
Getter for the Epetra_Import.
Definition: MapEpetra.cpp:394
double Real
Generic real data.
Definition: LifeV.hpp:175