LifeV
lifev/electrophysiology/testsuite/test_0DHodgkinHuxley/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 minimal model
30 
31  @date 01 - 2013
32  @author Simone Rossi <simone.rossi@epfl.ch>
33 
34  @contributor
35  @mantainer Simone Rossi <simone.rossi@epfl.ch>
36  */
37 
38 // Tell the compiler to ignore specific kind of warnings:
39 #pragma GCC diagnostic ignored "-Wunused-variable"
40 #pragma GCC diagnostic ignored "-Wunused-parameter"
41 
42 #include <Epetra_ConfigDefs.h>
43 #ifdef EPETRA_MPI
44 #include <mpi.h>
45 #include <Epetra_MpiComm.h>
46 #else
47 #include <Epetra_SerialComm.h>
48 #endif
49 
50 //Tell the compiler to restore the warning previously silented
51 #pragma GCC diagnostic warning "-Wunused-variable"
52 #pragma GCC diagnostic warning "-Wunused-parameter"
53 
54 
55 #include <fstream>
56 
57 #include <lifev/electrophysiology/solver/IonicModels/IonicHodgkinHuxley.hpp>
58 #include <lifev/core/LifeV.hpp>
59 
60 using namespace LifeV;
61 
62 #define SolutionTestNorm 18208.651394522181363
63 
64 
65 Int main ( Int argc, char** argv )
66 {
67  //********************************************//
68  // Import parameters from an xml list. Use //
69  // Teuchos to create a list from a given file //
70  // in the execution directory. //
71  //********************************************//
72  std::cout << "Importing parameters list...";
73  Teuchos::ParameterList parameterList = * ( Teuchos::getParametersFromXmlFile ( "HodgkinHuxleyParameters.xml" ) );
74  std::cout << " Done!" << std::endl;
75 
76  //********************************************//
77  // Creates a new model object representing the//
78  // model from Hodgkin Huxley model. //
79  //********************************************//
80  std::cout << "Building Constructor for Hodgkin Huxley Model with parameters ... ";
81  IonicHodgkinHuxley ionicModel (parameterList);
82  std::cout << " Done!" << std::endl;
83 
84  //********************************************//
85  // Show the parameters of the model as well as//
86  // other informations about the object. //
87  //********************************************//
88  ionicModel.showMe();
89 
90  //********************************************//
91  // Initialize the solution with the default //
92  // values //
93  //********************************************//
94  std::cout << "Initializing solution vector...";
95  std::vector<Real> states (ionicModel.Size(), 0);
96  ionicModel.initialize (states);
97  std::cout << " Done!" << std::endl;
98 
99  //********************************************//
100  // Initialize the rhs to 0. The rhs is the //
101  // vector containing the numerical values of //
102  // the time derivatives of the state //
103  // variables, that is, the right hand side of //
104  // the differential equation. //
105  //********************************************//
106  std::cout << "Initializing rhs..." ;
107  std::vector<Real> rhs (ionicModel.Size(), 0);
108  std::cout << " Done! " << std::endl;
109 
110  //********************************************//
111  // The model needs as external informations //
112  // the contraction velocity and the Calcium //
113  // concentration. //
114  //********************************************//
115  Real Iapp (0);
116 
117  //********************************************//
118  // Simulation starts on t=0 and ends on t=TF. //
119  // The timestep is given by dt //
120  //********************************************//
121  Real TF (40);
122  Real dt (0.005);
123  int useRushLarsen (1);
124 
125  //********************************************//
126  // Open the file "output.txt" to save the //
127  // solution. //
128  //********************************************//
129  std::string filename = "output.txt";
130  std::ofstream output ("output.txt");
131 
132  std::cout << "Potential: " << states.at (0) << std::endl;
133 
134  //********************************************//
135  // We record the norm of the solution to //
136  // check the failure of the test //
137  //********************************************//
138  Real SolutionNorm = states[0];
139 
140  //********************************************//
141  // Time loop starts. //
142  //********************************************//
143  std::cout << "Time loop starts...\n";
144  for ( Real t = 0; t < TF; )
145  {
146  //********************************************//
147  // Compute the applied current. This is a //
148  // simple switch. //
149  //********************************************//
150  if ( t > 20.5 && t < 21 )
151  {
152  Iapp = 380.0;
153  }
154  else
155  {
156  Iapp = 0;
157  }
158  std::cout << "\r " << t << " ms. " << std::flush;
159 
160  //********************************************//
161  // Compute the rhs using the model equations //
162  //********************************************//
163  ionicModel.setAppliedCurrent (Iapp);
164  if (useRushLarsen)
165  {
166  ionicModel.computeGatingVariablesWithRushLarsen (states, dt);
167  Real RHS = ionicModel.computeLocalPotentialRhs ( states);
168  ionicModel.addAppliedCurrent (RHS);
169  //********************************************//
170  // Use forward Euler method to advance the //
171  // solution in time. //
172  //********************************************//
173  states[0] = states[0] + dt * (RHS);
174  }
175  else
176  {
177  ionicModel.computeRhs ( states, rhs);
178  rhs[0] += Iapp;
179 
180  //********************************************//
181  // Use forward Euler method to advance the //
182  // solution in time. //
183  //********************************************//
184  for ( int j (0); j < ionicModel.Size(); j++)
185  {
186  states.at (j) = states.at (j) + dt * rhs.at (j);
187  }
188  }
189 
190 
191  //********************************************//
192  // Writes solution on file. //
193  //********************************************//
194  output << t << ", ";
195  for ( int j (0); j < ionicModel.Size() - 1; j++)
196  {
197  output << states.at (j) << ", ";
198  }
199  output << states.at ( ionicModel.Size() - 1 ) << "\n";
200 
201  //********************************************//
202  // Update the time. //
203  //********************************************//
204  t = t + dt;
205 
206  //********************************************//
207  // Update the norm of the solution to check //
208  // test failure //
209  //********************************************//
210  SolutionNorm += states[0];
211  }
212  std::cout << "\n...Time loop ends.\n";
213  std::cout << "Solution written on file: " << filename << "\n";
214 
215  //********************************************//
216  // Close exported file. //
217  //********************************************//
218  output.close();
219 
220  //********************************************//
221  // Check if the test failed //
222  //********************************************//
223  Real returnValue;
224 
225  Real err = std::abs (SolutionTestNorm - SolutionNorm) / std::abs (SolutionTestNorm);
226  std::cout << std::setprecision (20) << "\nError: " << err << "\nSolution norm: " << SolutionNorm << "\n";
227  if ( err > 1e-12 )
228  {
229  std::cout << "\nTest Failed!\n";
230  returnValue = EXIT_FAILURE; // Norm of solution did not match
231  }
232  else
233  {
234  returnValue = EXIT_SUCCESS;
235  }
236  return ( returnValue );
237 }
238 
239 
240 #undef SolutionTestNorm
int32_type Int
Generic integer data.
Definition: LifeV.hpp:188
void updateInverseJacobian(const UInt &iQuadPt)
int main(int argc, char **argv)
Definition: dummy.cpp:5
double Real
Generic real data.
Definition: LifeV.hpp:175