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