LifeV
lifev/electrophysiology/testsuite/test_0DLuoRudyIModel/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 Luo Rudy Phase I model
30 
31  @date 07 - 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/IonicLuoRudyI.hpp>
58 #include <lifev/core/LifeV.hpp>
59 
60 #include <Teuchos_RCP.hpp>
61 #include <Teuchos_ParameterList.hpp>
62 #include "Teuchos_XMLParameterListHelpers.hpp"
63 
64 #define SolutionTestNorm -14029.524312899517099
65 
66 using namespace LifeV;
67 
68 Int main ( Int argc, char** argv )
69 {
70  //! Initializing Epetra communicator
71  MPI_Init (&argc, &argv);
72  Epetra_MpiComm Comm (MPI_COMM_WORLD);
73  if ( Comm.MyPID() == 0 )
74  {
75  std::cout << "% using MPI" << std::endl;
76  }
77 
78 
79  //********************************************//
80  // Creates a new model object representing the//
81  // model from Luo Rudy Phase I model 1991. //
82  //********************************************//
83  std::cout << "Building Constructor for Luo Rudy I Model with parameters ... ";
84  IonicLuoRudyI ionicModel;
85  std::cout << " Done!" << std::endl;
86 
87 
88  //********************************************//
89  // Show the parameters of the model as well as//
90  // other informations about the object. //
91  //********************************************//
92  ionicModel.showMe();
93 
94 
95  //********************************************//
96  // Initialize the solution to default. //
97  //********************************************//
98  std::cout << "Initializing solution vector...";
99  std::vector<Real> states (ionicModel.restingConditions() );
100  std::cout << " Done!" << std::endl;
101 
102 
103  //********************************************//
104  // Initialize the rhs to 0. The rhs is the //
105  // vector containing the numerical values of //
106  // the time derivatives of the state //
107  // variables, that is, the right hand side of //
108  // the differential equation. //
109  //********************************************//
110  std::cout << "Initializing rhs..." ;
111  std::vector<Real> rhs (ionicModel.Size(), 0);
112  std::cout << " Done! " << std::endl;
113 
114 
115 
116 
117  //********************************************//
118  // The model needs as external informations //
119  // the contraction velocity and the Calcium //
120  // concentration. //
121  //********************************************//
122  Real Iapp (0.0);
123 
124 
125  //********************************************//
126  // Simulation starts on t=0 and ends on t=TF. //
127  // The timestep is given by dt //
128  //********************************************//
129  Real TF (500);
130  Real dt (.25);
131 
132  int savestep ( 1. / dt );
133  int iter (0);
134 
135  //********************************************//
136  // We record the norm of the solution to //
137  // check the failure of the test //
138  //********************************************//
139  Real SolutionNorm = states[0];
140 
141  //********************************************//
142  // Open the file "output.txt" to save the //
143  // solution. //
144  //********************************************//
145  std::string filename = "output.txt";
146  std::ofstream output ("output.txt");
147 
148  std::cout << "Potential: " << states[0] << std::endl;
149  //********************************************//
150  // Time loop starts. //
151  //********************************************//
152  std::cout << "Time loop starts...\n";
153  for ( Real t = 0; t < TF; )
154  {
155  iter++;
156  //********************************************//
157  // Compute the applied current. This is a //
158  // simple switch. //
159  //********************************************//
160  if ( t > 3 && t < 5 )
161  {
162  Iapp = 20.;
163  }
164  else
165  {
166  Iapp = 0;
167  }
168  std::cout << "\r " << t << " ms. " << std::flush;
169 
170  //********************************************//
171  // Compute the rhs using the model equations //
172  //********************************************//
173  ionicModel.setAppliedCurrent (Iapp);
174  ionicModel.computeRhs ( states, rhs);
175  ionicModel.addAppliedCurrent (rhs);
176 
177  //********************************************//
178  // Use forward Euler method to advance the //
179  // solution in time. //
180  //********************************************//
181 
182  states[0] = states[0] + dt * rhs[0];
183  ionicModel.computeGatingVariablesWithRushLarsen ( states, dt);
184 
185  int offset = 1 + ionicModel.numberOfGatingVariables();
186  for ( int j (0); j < ( ionicModel.Size() - offset ); j++)
187  {
188  states[j + offset] = states[j + offset] + dt * rhs[j + offset];
189  }
190 
191  //********************************************//
192  // Writes solution on file. //
193  //********************************************//
194  if ( iter % savestep == 0 )
195  {
196  output << t << ", ";
197  for ( int j (0); j < ionicModel.Size() - 1; j++)
198  {
199  output << states[j] << ", ";
200  }
201  output << states.at ( ionicModel.Size() - 1 ) << "\n";
202 
203  //********************************************//
204  // Update the norm of the solution to check //
205  // test failure //
206  //********************************************//
207  SolutionNorm = SolutionNorm + states[0];
208  }
209  //********************************************//
210  // Update the time. //
211  //********************************************//
212  t = t + dt;
213 
214 
215  }
216  std::cout << "\n...Time loop ends.\n";
217  std::cout << "Solution written on file: " << filename << "\n";
218  //********************************************//
219  // Close exported file. //
220  //********************************************//
221  output.close();
222 
223  //! Finalizing Epetra communicator
224  MPI_Finalize();
225  Real returnValue;
226 
227  Real err = std::abs (SolutionTestNorm - SolutionNorm) / std::abs (SolutionTestNorm);
228  std::cout << std::setprecision (20) << "\nError: " << err << "\nSolution norm: " << SolutionNorm << "\n";
229  if ( err > 1e-12 )
230  {
231  std::cout << "\nTest Failed!\n";
232  returnValue = EXIT_FAILURE; // Norm of solution did not match
233  }
234  else
235  {
236  returnValue = EXIT_SUCCESS;
237  }
238  return ( returnValue );
239 }
240 
241 #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