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/IonicHodgkinHuxley.hpp> 58 #include <lifev/core/LifeV.hpp> 60 using namespace LifeV;
62 #define SolutionTestNorm 18208.651394522181363
72 std::cout <<
"Importing parameters list...";
73 Teuchos::ParameterList parameterList = * ( Teuchos::getParametersFromXmlFile (
"HodgkinHuxleyParameters.xml" ) );
74 std::cout <<
" Done!" << std::endl;
80 std::cout <<
"Building Constructor for Hodgkin Huxley Model with parameters ... ";
81 IonicHodgkinHuxley ionicModel (parameterList);
82 std::cout <<
" Done!" << std::endl;
94 std::cout <<
"Initializing solution vector...";
95 std::vector<Real> states (ionicModel.Size(), 0);
96 ionicModel.initialize (states);
97 std::cout <<
" Done!" << std::endl;
106 std::cout <<
"Initializing rhs..." ;
107 std::vector<Real> rhs (ionicModel.Size(), 0);
108 std::cout <<
" Done! " << std::endl;
123 int useRushLarsen (1);
129 std::string filename =
"output.txt";
130 std::ofstream output (
"output.txt");
132 std::cout <<
"Potential: " << states.at (0) << std::endl;
138 Real SolutionNorm = states[0];
143 std::cout <<
"Time loop starts...\n";
144 for (
Real t = 0; t < TF; )
150 if ( t > 20.5 && t < 21 )
158 std::cout <<
"\r " << t <<
" ms. " << std::flush;
163 ionicModel.setAppliedCurrent (Iapp);
166 ionicModel.computeGatingVariablesWithRushLarsen (states, dt);
167 Real RHS = ionicModel.computeLocalPotentialRhs ( states);
168 ionicModel.addAppliedCurrent (RHS);
173 states[0] = states[0] + dt * (RHS);
177 ionicModel.computeRhs ( states, rhs);
184 for (
int j (0); j < ionicModel.Size(); j++)
186 states.at (j) = states.at (j) + dt * rhs.at (j);
195 for (
int j (0); j < ionicModel.Size() - 1; j++)
197 output << states.at (j) <<
", ";
199 output << states.at ( ionicModel.Size() - 1 ) <<
"\n";
210 SolutionNorm += states[0];
212 std::cout <<
"\n...Time loop ends.\n";
213 std::cout <<
"Solution written on file: " << filename <<
"\n";
226 std::cout << std::setprecision (20) <<
"\nError: " << err <<
"\nSolution norm: " << SolutionNorm <<
"\n";
229 std::cout <<
"\nTest Failed!\n";
230 returnValue = EXIT_FAILURE;
234 returnValue = EXIT_SUCCESS;
236 return ( returnValue );
240 #undef SolutionTestNorm
int32_type Int
Generic integer data.
void updateInverseJacobian(const UInt &iQuadPt)
int main(int argc, char **argv)
double Real
Generic real data.