LifeV
lifev/electrophysiology/testsuite/test_0DAlievPanfilovModel/main.cpp
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
29  @brief 0D test with the Aliev - Panfilov model 1996. Parameters are from Nash - Panfilov 2004.
30 
31  In this test we solve the Aliev - Panfilov ionic model alone as a set of ODEs.
32 
33  @date 01 - 2013
34  @author Simone Rossi <simone.rossi@epfl.ch>
35 
36  @contributor
37  @mantainer Simone Rossi <simone.rossi@epfl.ch>
38  */
39 
40 // To save the solution
41 #include <fstream>
42 
43 // To solve the model
44 #include <lifev/electrophysiology/solver/IonicModels/IonicAlievPanfilov.hpp>
45 
46 // To use LifeV
47 #include <lifev/core/LifeV.hpp>
48 
49 using namespace LifeV;
50 
51 //This is the norm of the precomputed solution
52 //we check the test against this value
53 #define SolutionTestNorm 3157.1639117380209427
54 
56 {
57  //********************************************//
58  // Import parameters from an xml list. Use //
59  // Teuchos to create a list from a given file //
60  // in the execution directory. //
61  //********************************************//
62 
63  std::cout << "Importing parameters list...";
64  Teuchos::ParameterList APParameterList = * ( Teuchos::getParametersFromXmlFile ( "AlievPanfilovParameters.xml" ) );
65  std::cout << " Done!" << std::endl;
66 
67 
68  //********************************************//
69  // Creates a new model object representing the//
70  // model from Aliev and Panfilov 1996. The //
71  // model input are the parameters. Pass the //
72  // parameter list in the constructor //
73  //********************************************//
74  std::cout << "Building Constructor for AlievPanfilov Model with parameters ... ";
75  IonicAlievPanfilov model ( APParameterList );
76  std::cout << " Done!" << std::endl;
77 
78 
79  //********************************************//
80  // Show the parameters of the model as well as//
81  // other informations about the object. //
82  //********************************************//
83  model.showMe();
84 
85 
86  //********************************************//
87  // Initialize the solution to 0. The model //
88  // consist of three state variables. //
89  // The method model.Size() //
90  // returns the number of state variables of //
91  // the model. rStates is the reference to the //
92  // the vector states //
93  //********************************************//
94  std::cout << "Initializing solution vector...";
95  std::vector<Real> unknowns (model.Size(), 0);
96  std::cout << " Done!" << std::endl;
97 
98  //********************************************//
99  // Initialize the rhs to 0. The rhs is the //
100  // vector containing the numerical values of //
101  // the time derivatives of the state //
102  // variables, that is, the right hand side of //
103  // the differential equation. //
104  //********************************************//
105  std::cout << "Initializing rhs..." ;
106  std::vector<Real> rhs (model.Size(), 0);
107  std::cout << " Done! " << std::endl;
108 
109  //********************************************//
110  // The model needs as external informations //
111  // the applied current. //
112  //********************************************//
113  Real Iapp (0.0);
114 
115  //********************************************//
116  // Simulation starts on t=0 and ends on t=TF. //
117  // The timestep is given by dt //
118  //********************************************//
119  Real TF (60);
120  Real dt (0.01);
121 
122  //********************************************//
123  // Open the file "output.txt" to save the //
124  // solution. //
125  //********************************************//
126  std::string filename = "output.txt";
127  std::ofstream output ("output.txt");
128 
129  //********************************************//
130  // We record the norm of the solution to //
131  // check the failure of the test //
132  //********************************************//
133  Real SolutionNorm = unknowns[0];
134 
135  //********************************************//
136  // Time loop starts. //
137  //********************************************//
138  std::cout << "Time loop starts...\n";
139  for ( Real t = 0; t < TF; )
140  {
141 
142  //********************************************//
143  // Compute the applied current. This is a //
144  // simple switch. //
145  //********************************************//
146  if ( t > 0.1 && t < 0.2 )
147  {
148  Iapp = 2.0;
149  }
150  else
151  {
152  Iapp = 0;
153  }
154 
155  //********************************************//
156  // Set the applied current to the model. //
157  //********************************************//
158  model.setAppliedCurrent (Iapp);
159  std::cout << "\r " << t << " ms. " << std::flush;
160 
161  //********************************************//
162  // Compute the rhs using the model equations //
163  //********************************************//
164  model.computeRhs ( unknowns, rhs);
165  model.addAppliedCurrent (rhs);
166 
167  //********************************************//
168  // Use forward Euler method to advance the //
169  // solution in time. //
170  //********************************************//
171  unknowns.at (0) = unknowns.at (0) + dt * rhs.at (0);
172  unknowns.at (1) = unknowns.at (1) + dt * rhs.at (1);
173 
174  //********************************************//
175  // Writes solution on file. //
176  //********************************************//
177  output << t << ", " << unknowns.at (0) << ", " << unknowns.at (1) << "\n";
178 
179  //********************************************//
180  // Update the time. //
181  //********************************************//
182  t = t + dt;
183 
184  //********************************************//
185  // Update the norm of the solution to check //
186  // test failure //
187  //********************************************//
188  SolutionNorm += unknowns[0];
189  }
190  std::cout << "\n...Time loop ends.\n";
191  std::cout << "Solution written on file: " << filename << "\n";
192  //********************************************//
193  // Close exported file. //
194  //********************************************//
195  output.close();
196 
197  Real returnValue;
198  Real err = std::abs (SolutionTestNorm - SolutionNorm) / std::abs (SolutionTestNorm);
199  std::cout << std::setprecision (20) << "\nError: " << err << "\nSolution norm: " << SolutionNorm << "\n";
200  if ( err > 1e-12 )
201  {
202  returnValue = EXIT_FAILURE; // Norm of solution did not match
203  }
204  else
205  {
206  returnValue = EXIT_SUCCESS;
207  }
208  return ( returnValue );
209 }
210 
211 #undef SolutionTestNorm
int32_type Int
Generic integer data.
Definition: LifeV.hpp:188
void updateInverseJacobian(const UInt &iQuadPt)
double Real
Generic real data.
Definition: LifeV.hpp:175