39 #pragma GCC diagnostic ignored "-Wunused-variable" 40 #pragma GCC diagnostic ignored "-Wunused-parameter" 42 #include <Epetra_ConfigDefs.h> 45 #include <Epetra_MpiComm.h> 47 #include <Epetra_SerialComm.h> 51 #pragma GCC diagnostic warning "-Wunused-variable" 52 #pragma GCC diagnostic warning "-Wunused-parameter" 57 #include <lifev/electrophysiology/solver/IonicModels/IonicFox.hpp> 58 #include <lifev/core/LifeV.hpp> 60 #include <Teuchos_RCP.hpp> 61 #include <Teuchos_ParameterList.hpp> 62 #include "Teuchos_XMLParameterListHelpers.hpp" 64 using namespace LifeV;
68 #define SolutionTestNorm -593906.94845120306127
73 MPI_Init (&argc, &argv);
74 Epetra_MpiComm Comm (MPI_COMM_WORLD);
75 if ( Comm.MyPID() == 0 )
77 std::cout <<
"% using MPI" << std::endl;
86 std::cout <<
"Importing parameters list...";
87 Teuchos::ParameterList parameterList = * ( Teuchos::getParametersFromXmlFile (
"FoxParameters.xml" ) );
88 std::cout <<
" Done!" << std::endl;
96 std::cout <<
"Building Constructor for Fox Model with parameters ... ";
98 std::cout <<
" Done!" << std::endl;
110 std::cout <<
"Initializing solution vector...";
111 std::vector<Real> unknowns (model.Size(), 0 );
112 model.initialize (unknowns);
113 std::cout <<
" Done!" << std::endl;
122 std::cout <<
"Initializing rhs..." ;
123 std::vector<Real> rhs (model.Size() + 1, 0);
124 std::cout <<
" Done! " << std::endl;
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 );
146 std::string filename =
"output.txt";
147 std::ofstream output (
"output.txt");
153 std::cout <<
"Time loop starts...\n";
157 int savedt ( parameterList.get (
"savedt", 1.0) / dt );
163 Real SolutionNorm = unknowns[0];
165 for (
Real t = 0; t < TF; )
172 if ( t >= timeSt && t <= timeSt + 1.0 )
175 if ( t >= timeSt + 1.0 - dt && t <= timeSt + 1.0 )
177 timeSt = timeSt + stInt;
185 std::cout <<
"\r " << t <<
" ms. " << std::flush;
190 model.setAppliedCurrent (Iapp);
191 model.computeRhs ( unknowns, rhs );
192 model.addAppliedCurrent (rhs);
200 if ( iter % savedt == 0)
206 SolutionNorm += unknowns[0];
207 output << t <<
", " << unknowns.at (0) <<
", " << unknowns.at (10) <<
", " << unknowns.at (11)
208 <<
", " << rhs.at (13) <<
"\n";
217 for (
int j (0); j <= 12; ++j)
219 unknowns.at (j) = unknowns.at (j) + dt * rhs.at (j);
228 std::cout <<
"\n...Time loop ends.\n";
229 std::cout <<
"Solution written on file: " << filename <<
"\n";
241 std::cout << std::setprecision (20) <<
"\nError: " << err <<
"\nSolution norm: " << SolutionNorm <<
"\n";
245 std::cout <<
"\nTest Failed!\n";
246 returnValue = EXIT_FAILURE;
250 returnValue = EXIT_SUCCESS;
252 return ( returnValue );
256 #undef SolutionTestNorm void showMe()
Display information about the model.
int32_type Int
Generic integer data.
void updateInverseJacobian(const UInt &iQuadPt)
int main(int argc, char **argv)
IonicModel - This class implements an ionic model.
double Real
Generic real data.