LifeV
lifev/electrophysiology/testsuite/test_0DNoblePurkinje/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 09 - 2013
32  @author Simone Palamara <palamara.simone@gmail.com>
33 
34  @contributor Simone Rossi <simone.rossi@epfl.ch>
35  */
36 
37 // Tell the compiler to ignore specific kind of warnings:
38 #pragma GCC diagnostic ignored "-Wunused-variable"
39 #pragma GCC diagnostic ignored "-Wunused-parameter"
40 
41 #include <Epetra_ConfigDefs.h>
42 #ifdef EPETRA_MPI
43 #include <mpi.h>
44 #include <Epetra_MpiComm.h>
45 #else
46 #include <Epetra_SerialComm.h>
47 #endif
48 
49 //Tell the compiler to restore the warning previously silented
50 #pragma GCC diagnostic warning "-Wunused-variable"
51 #pragma GCC diagnostic warning "-Wunused-parameter"
52 
53 
54 #include <fstream>
55 
56 #include <lifev/electrophysiology/solver/IonicModels/IonicNoblePurkinje.hpp>
57 #include <lifev/core/LifeV.hpp>
58 
59 #include <Teuchos_RCP.hpp>
60 #include <Teuchos_ParameterList.hpp>
61 #include "Teuchos_XMLParameterListHelpers.hpp"
62 
63 using namespace LifeV;
64 
65 #define SolutionTestNorm -932355.10562017653137
66 
67 Int main ( Int argc, char** argv )
68 {
69  //! Initializing Epetra communicator
70  MPI_Init (&argc, &argv);
71  Epetra_MpiComm Comm (MPI_COMM_WORLD);
72  if ( Comm.MyPID() == 0 )
73  {
74  std::cout << "% using MPI" << std::endl;
75  }
76 
77  //********************************************//
78  // Creates a new model object representing the//
79  // model from Negroni and Lascano 1996. The //
80  // model input are the parameters. Pass the //
81  // parameter list in the constructor //
82  //********************************************//
83  std::cout << "Building Constructor for Noble Model with parameters ... ";
84  IonicNoblePurkinje 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  // Initialize the solution with the default //
96  // values //
97  //********************************************//
98  std::cout << "Initializing solution vector...";
99  std::vector<Real> states (ionicModel.restingConditions() );
100  ionicModel.initialize (states);
101  std::cout << " Done!" << std::endl;
102 
103 
104  //********************************************//
105  // Initialize the rhs to 0. The rhs is the //
106  // vector containing the numerical values of //
107  // the time derivatives of the state //
108  // variables, that is, the right hand side of //
109  // the differential equation. //
110  //********************************************//
111  std::cout << "Initializing rhs..." ;
112  std::vector<Real> rhs (ionicModel.Size(), 0);
113  std::cout << " Done! " << std::endl;
114 
115 
116 
117 
118  //********************************************//
119  // The model needs as external informations //
120  // the contraction velocity and the Calcium //
121  // concentration. //
122  //********************************************//
123  Real Iapp (0);
124 
125 
126  //********************************************//
127  // Simulation starts on t=0 and ends on t=TF. //
128  // The timestep is given by dt //
129  //********************************************//
130  Real TF (2000);
131  Real dt (0.1);
132  int useRushLarsen (1);
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  // Open the file "output.txt" to save the //
142  // solution. //
143  //********************************************//
144  std::string filename = "output.txt";
145  std::ofstream output ("output.txt");
146 
147  std::cout << "Potential: " << states.at (0) << std::endl;
148  //********************************************//
149  // Time loop starts. //
150  //********************************************//
151  std::cout << "Time loop starts...\n";
152  for ( Real t = 0; t < TF; )
153  {
154 
155  //********************************************//
156  // Compute Calcium concentration. Here it is //
157  // given as a function of time. //
158  //********************************************//
159  if ( t > 50.5 && t < 52 )
160  {
161  Iapp = 10.0;
162  }
163  else
164  {
165  Iapp = 0;
166  }
167  // Iapp = 0;
168  std::cout << "\r " << t << " ms. " << std::flush;
169 
170  //********************************************//
171  // Compute the rhs using the model equations //
172  //********************************************//
173  ionicModel.setAppliedCurrent (Iapp);
174  if (useRushLarsen)
175  {
176  ionicModel.computeGatingVariablesWithRushLarsen (states, dt);
177  Real RHS = ionicModel.computeLocalPotentialRhs ( states);
178  RHS += Iapp;
179  //********************************************//
180  // Use forward Euler method to advance the //
181  // solution in time. //
182  //********************************************//
183 
184  for ( int j (0); j < 1; j++)
185  {
186  states.at (j) = states.at (j) + dt * (RHS);
187  }
188  }
189  else
190  {
191  ionicModel.computeRhs ( states, rhs);
192  rhs[0] += Iapp;
193  //********************************************//
194  // Use forward Euler method to advance the //
195  // solution in time. //
196  //********************************************//
197 
198  for ( int j (0); j < ionicModel.Size(); j++)
199  {
200  states.at (j) = states.at (j) + dt * rhs.at (j);
201  }
202  }
203 
204  //********************************************//
205  // Writes solution on file. //
206  //********************************************//
207  output << t << ", ";
208  for ( int j (0); j < ionicModel.Size() - 1; j++)
209  {
210  output << states.at (j) << ", ";
211  }
212  output << states.at ( ionicModel.Size() - 1 ) << "\n";
213 
214  //********************************************//
215  // Update the time. //
216  //********************************************//
217  t = t + dt;
218 
219  //********************************************//
220  // Update the norm of the solution to check //
221  // test failure //
222  //********************************************//
223  SolutionNorm += states[0];
224  }
225  std::cout << "\n...Time loop ends.\n";
226  std::cout << "Solution written on file: " << filename << "\n";
227  //********************************************//
228  // Close exported file. //
229  //********************************************//
230  output.close();
231 
232  //! Finalizing Epetra communicator
233  MPI_Finalize();
234  Real returnValue;
235 
236  Real err = std::abs (SolutionTestNorm - SolutionNorm) / std::abs (SolutionTestNorm);
237  std::cout << std::setprecision (20) << "\nError: " << err << "\nSolution norm: " << SolutionNorm << "\n";
238  if ( err > 1e-12 )
239  {
240  std::cout << "\nTest Failed!\n";
241  returnValue = EXIT_FAILURE; // Norm of solution did not match
242  }
243  else
244  {
245  returnValue = EXIT_SUCCESS;
246  }
247  return ( returnValue );
248 }
249 
250 #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