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