LifeV
lifev/electrophysiology/testsuite/test_0DFoxModel/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 Fox model
30 
31  @date 04 - 2013
32  @author Marie Dupraz <dupraz.marie@gmail.com>
33 
34  @contributor
35  @mantainer Marie Dupraz <dupraz.marie@gmail.com>
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/IonicFox.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 using namespace LifeV;
65 
66 //This is the norm of the precomputed solution
67 //we check the test against this value
68 #define SolutionTestNorm -593906.94845120306127
69 
70 Int main ( Int argc, char** argv )
71 {
72  //! Initializing Epetra communicator
73  MPI_Init (&argc, &argv);
74  Epetra_MpiComm Comm (MPI_COMM_WORLD);
75  if ( Comm.MyPID() == 0 )
76  {
77  std::cout << "% using MPI" << std::endl;
78  }
79 
80  //********************************************//
81  // Import parameters from an xml list. Use //
82  // Teuchos to create a list from a given file //
83  // in the execution directory. //
84  //********************************************//
85 
86  std::cout << "Importing parameters list...";
87  Teuchos::ParameterList parameterList = * ( Teuchos::getParametersFromXmlFile ( "FoxParameters.xml" ) );
88  std::cout << " Done!" << std::endl;
89 
90  //********************************************//
91  // Creates a new model object representing the//
92  // model from Fox 2002. The //
93  // model input are the parameters. Pass the //
94  // parameter list in the constructor //
95  //********************************************//
96  std::cout << "Building Constructor for Fox Model with parameters ... ";
97  IonicFox model;
98  std::cout << " Done!" << std::endl;
99 
100  //********************************************//
101  // Show the parameters of the model as well as//
102  // other informations about the object. //
103  //********************************************//
104  model.showMe();
105 
106  //********************************************//
107  // Initialize the solution with the default //
108  // values //
109  //********************************************//
110  std::cout << "Initializing solution vector...";
111  std::vector<Real> unknowns (model.Size(), 0 );
112  model.initialize (unknowns);
113  std::cout << " Done!" << std::endl;
114 
115  //********************************************//
116  // Initialize the rhs to 0. The rhs is the //
117  // vector containing the numerical values of //
118  // the time derivatives of the state //
119  // variables, that is, the right hand side of //
120  // the differential equation. //
121  //********************************************//
122  std::cout << "Initializing rhs..." ;
123  std::vector<Real> rhs (model.Size() + 1, 0);
124  std::cout << " Done! " << std::endl;
125 
126  //********************************************//
127  // The model needs as external informations //
128  // the contraction velocity and the Calcium //
129  // concentration. //
130  //********************************************//
131  Real Iapp (0.0);
132 
133  //********************************************//
134  // Simulation starts on t=0 and ends on t=TF. //
135  // The timestep is given by dt //
136  //********************************************//
137  Real TF = parameterList.get ( "endTime", 5.0 );
138  Real dt = parameterList.get ( "timeStep", 0.005 );
139  Real timeSt = parameterList.get ( "stimuliTime", 1.0 );
140  Real stInt = parameterList.get ( "stimuliInterval", 1000.0 );
141 
142  //********************************************//
143  // Open the file "output.txt" to save the //
144  // solution. //
145  //********************************************//
146  std::string filename = "output.txt";
147  std::ofstream output ("output.txt");
148 
149 
150  //********************************************//
151  // Time loop starts. //
152  //********************************************//
153  std::cout << "Time loop starts...\n";
154 
155 
156  int iter (0);
157  int savedt ( parameterList.get ( "savedt", 1.0) / dt );
158 
159  //********************************************//
160  // We record the norm of the solution to //
161  // check the failure of the test //
162  //********************************************//
163  Real SolutionNorm = unknowns[0];
164 
165  for ( Real t = 0; t < TF; )
166  {
167 
168  //********************************************//
169  // Compute the applied current. This is a //
170  // simple switch. //
171  //********************************************//
172  if ( t >= timeSt && t <= timeSt + 1.0 )
173  {
174  Iapp = 80.; //0.516289;
175  if ( t >= timeSt + 1.0 - dt && t <= timeSt + 1.0 )
176  {
177  timeSt = timeSt + stInt;
178  }
179  }
180  else
181  {
182  Iapp = 0.;
183  }
184 
185  std::cout << "\r " << t << " ms. " << std::flush;
186 
187  //********************************************//
188  // Compute the rhs using the model equations //
189  //********************************************//
190  model.setAppliedCurrent (Iapp);
191  model.computeRhs ( unknowns, rhs );
192  model.addAppliedCurrent (rhs);
193 
194 
195  //********************************************//
196  // Writes solution on file. //
197  //********************************************//
198 
199 
200  if ( iter % savedt == 0)
201  {
202  //********************************************//
203  // Update the norm of the solution to check //
204  // test failure //
205  //********************************************//
206  SolutionNorm += unknowns[0];
207  output << t << ", " << unknowns.at (0) << ", " << unknowns.at (10) << ", " << unknowns.at (11)
208  << ", " << rhs.at (13) << "\n"; //
209 
210  }
211 
212  //********************************************//
213  // Use forward Euler method to advance the //
214  // solution in time. //
215  //********************************************//
216 
217  for (int j (0); j <= 12; ++j)
218  {
219  unknowns.at (j) = unknowns.at (j) + dt * rhs.at (j);
220  }
221 
222  //********************************************//
223  // Update the time. //
224  //********************************************//
225  t = t + dt;
226  }
227 
228  std::cout << "\n...Time loop ends.\n";
229  std::cout << "Solution written on file: " << filename << "\n";
230 
231  //********************************************//
232  // Close exported file. //
233  //********************************************//
234  output.close();
235 
236 
237  //! Finalizing Epetra communicator
238  MPI_Finalize();
239  Real returnValue;
240  Real err = std::abs (SolutionTestNorm - SolutionNorm) / std::abs (SolutionTestNorm);
241  std::cout << std::setprecision (20) << "\nError: " << err << "\nSolution norm: " << SolutionNorm << "\n";
242  if ( err > 1e-12 )
243 
244  {
245  std::cout << "\nTest Failed!\n";
246  returnValue = EXIT_FAILURE; // Norm of solution did not match
247  }
248  else
249  {
250  returnValue = EXIT_SUCCESS;
251  }
252  return ( returnValue );
253 }
254 
255 
256 #undef SolutionTestNorm
void showMe()
Display information about the model.
Definition: IonicFox.cpp:640
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
IonicModel - This class implements an ionic model.
Definition: IonicFox.hpp:61
double Real
Generic real data.
Definition: LifeV.hpp:175