void TestWithHeterogeneousCellModels() { HeartConfig::Instance()->SetSimulationDuration(1.0); //ms HeartConfig::Instance()->SetUseStateVariableInterpolation(true); HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.005, 0.01, 1.0); TetrahedralMesh<1,1> mesh; mesh.ConstructRegularSlabMesh(0.01, 1.0); HeterogeneousCellFactory cell_factory; MonodomainProblem<1> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); //// It's really very difficult to get hold of the correction assembler from here. //// The following, which requires 2 classes to be friends of this test compiles //// but fails asserts in the first line as bcc is not set up. //AbstractDynamicLinearPdeSolver<1,1,1>* p_solver = monodomain_problem.CreateSolver(); //MonodomainSolver<1,1>* p_mono_solver = dynamic_cast<MonodomainSolver<1,1>*>(p_solver); //MonodomainCorrectionTermAssembler<1,1>* p_assembler = p_mono_solver->mpMonodomainCorrectionTermAssembler; //TS_ASSERT_EQUALS(p_assembler->mElementsCanDoSvi.size(), 10u); // Therefore, we just test that calling Solve() runs (without the checking that // cell models are identical, this fails). monodomain_problem.Solve(); }
void TestCompareRealisticGeometry() throw(Exception) { HeartConfig::Reset(); HeartConfig::Instance()->SetSimulationDuration(50); //ms HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.01,0.1,0.1); double spatial_step = 0.05; HeartConfig::Instance()->SetMeshFileName("heart/test/data/UCSD_heart"); /* * Standard solve */ HeartConfig::Instance()->SetOutputDirectory("CompareRealisticGeometryStandard"); HeartConfig::Instance()->SetOutputFilenamePrefix("CompareRealisticGeometryStandard"); SmallBenchmarkStimulusHeartCellFactory<CellLuoRudy1991FromCellMLBackwardEuler> cell_factory; MonodomainProblem<3> monodomain_problem( &cell_factory ); monodomain_problem.Initialise(); monodomain_problem.Solve(); DistributedVector standard_solution = monodomain_problem.GetSolutionDistributedVector(); HeartEventHandler::Headings(); HeartEventHandler::Report(); /* * Mass lumping solve */ HeartEventHandler::Reset(); HeartConfig::Instance()->SetOutputDirectory("CompareRealisticGeometryMassLumping"); HeartConfig::Instance()->SetOutputFilenamePrefix("CompareRealisticGeometryMassLumping"); HeartConfig::Instance()->SetUseMassLumping(); MonodomainProblem<3> monodomain_problem_ml( &cell_factory ); monodomain_problem_ml.Initialise(); monodomain_problem_ml.Solve(); DistributedVector mass_lumping_solution = monodomain_problem_ml.GetSolutionDistributedVector(); HeartEventHandler::Headings(); HeartEventHandler::Report(); // The idea is to check that the error stays O(h) double tolerance = 100*spatial_step; for (DistributedVector::Iterator index = standard_solution.Begin(); index != standard_solution.End(); ++index) { TS_ASSERT_DELTA(standard_solution[index], mass_lumping_solution[index], tolerance); } }
// Solve on a 1D string of cells, 1mm long with a space step of 0.1mm. void TestMonodomainConstantStimulus() { // this parameters are a bit arbitrary, and chosen to get a good spread of voltages HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(1.75)); HeartConfig::Instance()->SetSimulationDuration(2); //ms HeartConfig::Instance()->SetMeshFileName("mesh/test/data/1D_0_to_1mm_10_elements"); HeartConfig::Instance()->SetOutputDirectory("MonoNeumannConst"); HeartConfig::Instance()->SetOutputFilenamePrefix("MonodomainLR91_1d"); ZeroStimulusCellFactory<CellLuoRudy1991FromCellML, 1> cell_factory; MonodomainProblem<1> monodomain_problem( &cell_factory ); monodomain_problem.Initialise(); HeartConfig::Instance()->SetSurfaceAreaToVolumeRatio(1*1.75/0.0005); // create boundary conditions container boost::shared_ptr<BoundaryConditionsContainer<1,1,1> > p_bcc(new BoundaryConditionsContainer<1,1,1>); ConstBoundaryCondition<1>* p_bc_stim = new ConstBoundaryCondition<1>(2*1.75/0.0005); // get mesh AbstractTetrahedralMesh<1,1> &mesh = monodomain_problem.rGetMesh(); // loop over boundary elements AbstractTetrahedralMesh<1, 1>::BoundaryElementIterator iter; iter = mesh.GetBoundaryElementIteratorBegin(); while (iter != mesh.GetBoundaryElementIteratorEnd()) { // if the element is on the left of the mesh, add a stimulus to the bcc if (((*iter)->GetNodeLocation(0))[0]==0.0) { p_bcc->AddNeumannBoundaryCondition(*iter, p_bc_stim); } iter++; } // pass the bcc to the monodomain problem monodomain_problem.SetBoundaryConditionsContainer(p_bcc); monodomain_problem.Solve(); // check some voltages ReplicatableVector voltage_replicated(monodomain_problem.GetSolution()); double atol=5e-3; TS_ASSERT_DELTA(voltage_replicated[1], 94.6426, atol); TS_ASSERT_DELTA(voltage_replicated[3], 49.7867, atol); TS_ASSERT_DELTA(voltage_replicated[5], 30.5954, atol); TS_ASSERT_DELTA(voltage_replicated[7], 21.6782, atol); TS_ASSERT_DELTA(voltage_replicated[9], -33.9983, atol); TS_ASSERT_DELTA(voltage_replicated[10], -52.2396, atol); }
void TestCoverage3d() { HeartConfig::Instance()->SetSimulationDuration(0.1); //ms HeartConfig::Instance()->SetUseStateVariableInterpolation(true); HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.005, 0.01, 0.1); TetrahedralMesh<3,3> mesh; mesh.ConstructRegularSlabMesh(0.02, 0.02, 0.02, 0.02); ZeroStimulusCellFactory<CellLuoRudy1991FromCellML,3> cell_factory; MonodomainProblem<3> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); }
//////////////////////////////////////////////////////////// // Compare Mono and Bidomain Simulations //////////////////////////////////////////////////////////// void TestCompareBidomainProblemWithMonodomain3D() { // the bidomain equations reduce to the monodomain equations // if sigma_e is infinite (equivalent to saying the extra_cellular // space is grounded. sigma_e is set to be very large here: HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(1.75, 1.75, 1.75)); HeartConfig::Instance()->SetExtracellularConductivities(Create_c_vector(17500, 17500, 17500)); HeartConfig::Instance()->SetSimulationDuration(1.0); //ms HeartConfig::Instance()->SetMeshFileName("mesh/test/data/3D_0_to_1mm_6000_elements"); HeartConfig::Instance()->SetOutputDirectory("Monodomain3d"); HeartConfig::Instance()->SetOutputFilenamePrefix("monodomain3d"); /////////////////////////////////////////////////////////////////// // monodomain /////////////////////////////////////////////////////////////////// PlaneStimulusCellFactory<CellLuoRudy1991FromCellML, 3> cell_factory(-600.0*1000); MonodomainProblem<3> monodomain_problem( &cell_factory ); monodomain_problem.Initialise(); monodomain_problem.Solve(); /////////////////////////////////////////////////////////////////// // bidomain /////////////////////////////////////////////////////////////////// HeartConfig::Instance()->SetOutputDirectory("Bidomain3d"); HeartConfig::Instance()->SetOutputFilenamePrefix("bidomain3d"); BidomainProblem<3> bidomain_problem( &cell_factory ); bidomain_problem.Initialise(); bidomain_problem.Solve(); /////////////////////////////////////////////////////////////////// // compare /////////////////////////////////////////////////////////////////// DistributedVector monodomain_voltage = monodomain_problem.GetSolutionDistributedVector(); DistributedVector bidomain_solution = bidomain_problem.GetSolutionDistributedVector(); DistributedVector::Stripe bidomain_voltage(bidomain_solution,0); DistributedVector::Stripe extracellular_potential(bidomain_solution,1); for (DistributedVector::Iterator index = bidomain_solution.Begin(); index != bidomain_solution.End(); ++index) { TS_ASSERT_DELTA(monodomain_voltage[index], bidomain_voltage[index], 0.5); TS_ASSERT_DELTA(extracellular_potential[index], 0, 1.0); } }
/* * This is the same as TestConductionVelocityConvergesFasterWithSvi1d with i=2, but solves in two parts. * If that test changes, check the hardcoded values here! */ void TestArchiving() { std::string archive_dir = "monodomain_svi_archive"; std::string archive_file = "monodomain_svi.arch"; std::string output_dir = "monodomain_svi_output"; { // Save HeartConfig::Instance()->SetSimulationDuration(0.1); //ms HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.01, 0.01, 0.01); HeartConfig::Instance()->SetOutputDirectory(output_dir); HeartConfig::Instance()->SetOutputFilenamePrefix("results"); HeartConfig::Instance()->SetUseStateVariableInterpolation(); DistributedTetrahedralMesh<1,1> mesh; mesh.ConstructRegularSlabMesh(0.02, 1.0); BlockCellFactory<1> cell_factory; MonodomainProblem<1> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); CardiacSimulationArchiver<MonodomainProblem<1> >::Save(monodomain_problem, archive_dir); } HeartConfig::Instance()->Reset(); { // Load HeartConfig::Instance()->SetUseStateVariableInterpolation(false); // Just in case... MonodomainProblem<1> *p_monodomain_problem = CardiacSimulationArchiver<MonodomainProblem<1> >::Load(archive_dir); HeartConfig::Instance()->SetSimulationDuration(4.0); //ms p_monodomain_problem->Solve(); ReplicatableVector final_voltage; final_voltage.ReplicatePetscVector(p_monodomain_problem->GetSolution()); double probe_voltage = final_voltage[15]; TS_ASSERT_DELTA(probe_voltage, 17.3131, 1e-3); delete p_monodomain_problem; } }
// Solve with a space step of 0.5mm. // // Note that this space step ought to be too big! // // NOTE: This test uses NON-PHYSIOLOGICAL parameters values (conductivities, // surface-area-to-volume ratio, capacitance, stimulus amplitude). Essentially, // the equations have been divided through by the surface-area-to-volume ratio. // (Historical reasons...) void TestMonodomainFailing() { if (PetscTools::GetNumProcs() > 3u) { TS_TRACE("This test is not suitable for more than 3 processes."); return; } HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(0.0005)); HeartConfig::Instance()->SetSimulationDuration(1); //ms // HeartConfig::Instance()->SetMeshFileName("mesh/test/data/1D_0_to_1_20_elements"); // HeartConfig::Instance()->SetSpaceDimension(1); // HeartConfig::Instance()->SetFibreLength(1.0, 1.0/20.0); HeartConfig::Instance()->SetSpaceDimension(3); HeartConfig::Instance()->SetSlabDimensions(0.2, 0.1, 0.1, 0.05); HeartConfig::Instance()->SetOutputDirectory("MonoFailing"); HeartConfig::Instance()->SetOutputFilenamePrefix("MonodomainLR91_SVI"); HeartConfig::Instance()->SetUseStateVariableInterpolation(true); PlaneStimulusCellFactory<CellLuoRudy1991FromCellML, 3> cell_factory; MonodomainProblem<3> monodomain_problem(&cell_factory); monodomain_problem.Initialise(); HeartConfig::Instance()->SetSurfaceAreaToVolumeRatio(1.0); // HeartConfig::Instance()->SetCapacitance(1.0); // the mesh is too coarse, and this simulation will result in cell gating // variables going out of range. An exception should be thrown in the // EvaluateYDerivatives() method of the cell model TS_ASSERT_THROWS_CONTAINS(monodomain_problem.Solve(), #ifndef NDEBUG // VerifyStateVariables is only called when debug is on "State variable fast_sodium_current_m_gate__m has gone out of range. " "Check numerical parameters, for example time and space stepsizes"); #else // This test hits a later assert(!isnan) if we do a ndebug build. "Assertion tripped: !std::isnan(i_ionic)"); #endif // NDEBUG }
void TestConductionVelocityInCrossFibreDirection2d() { ReplicatableVector final_voltage_ici; ReplicatableVector final_voltage_svi; ReplicatableVector final_voltage_svit; HeartConfig::Instance()->SetSimulationDuration(5.0); //ms //HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.0005, 0.01, 0.01); //See comment below HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.005, 0.01, 0.01); // much lower conductivity in cross-fibre direction - ICI will struggle HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(1.75, 0.17)); TetrahedralMesh<2,2> mesh; mesh.ConstructRegularSlabMesh(0.02 /*h*/, 0.5, 0.3); // ICI - nodal current interpolation - the default { HeartConfig::Instance()->SetOutputDirectory("MonodomainIci2d"); HeartConfig::Instance()->SetOutputFilenamePrefix("results"); HeartConfig::Instance()->SetUseStateVariableInterpolation(false); BlockCellFactory<2> cell_factory; MonodomainProblem<2> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); final_voltage_ici.ReplicatePetscVector(monodomain_problem.GetSolution()); } // SVI - state variable interpolation { HeartConfig::Instance()->SetOutputDirectory("MonodomainSvi2d"); HeartConfig::Instance()->SetOutputFilenamePrefix("results"); HeartConfig::Instance()->SetUseStateVariableInterpolation(); BlockCellFactory<2> cell_factory; MonodomainProblem<2> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); final_voltage_svi.ReplicatePetscVector(monodomain_problem.GetSolution()); } #ifdef CHASTE_CVODE ReplicatableVector final_voltage_svi_cvode; // SVI - state variable interpolation with CVODE cells { HeartConfig::Instance()->SetOutputDirectory("MonodomainSvi2dCvode"); HeartConfig::Instance()->SetOutputFilenamePrefix("results"); HeartConfig::Instance()->SetUseStateVariableInterpolation(); BlockCellFactoryCvode<2> cell_factory; MonodomainProblem<2> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); final_voltage_svi_cvode.ReplicatePetscVector(monodomain_problem.GetSolution()); } #endif //CHASTE_CVODE // SVIT - state variable interpolation on non-distributed tetrahedral mesh { HeartConfig::Instance()->SetOutputDirectory("MonodomainSviTet2d"); HeartConfig::Instance()->SetOutputFilenamePrefix("results"); HeartConfig::Instance()->SetUseStateVariableInterpolation(); BlockCellFactory<2> cell_factory; MonodomainProblem<2> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); final_voltage_svit.ReplicatePetscVector(monodomain_problem.GetSolution()); } // Visualised results with h=0.02 and h=0.01 - results looks sensible according to // paper: // 1. SVI h=0.01 and h=0.02 match more closely than ICI results - ie SVI converges faster // 2. CV in fibre direction faster for ICI (both values of h) // 3. CV in cross fibre direction: (i) h=0.01, faster for ICI; h=0.02 slower for ICI. // (Matches results in paper) double ici_20 = -17.1939; // These numbers are from a solve with ODE timestep of 0.0005, double svi_20 = -62.6336; // i.e. ten times smaller than that used in the test at present double ici_130 = 15.4282; // (for speed), hence large tolerances below. double svi_130 = 30.7389; // The tolerances are still nowhere near overlapping - i.e. ICI different to SVI // node 20 (for h=0.02) is on the x-axis (fibre direction), SVI CV is slower TS_ASSERT_DELTA(mesh.GetNode(20)->rGetLocation()[0], 0.4, 1e-9); TS_ASSERT_DELTA(mesh.GetNode(20)->rGetLocation()[1], 0.0, 1e-9); TS_ASSERT_DELTA(final_voltage_ici[20], ici_20, 8.0); // These tolerances show difference in parallel, TS_ASSERT_DELTA(final_voltage_svi[20], svi_20, 3.0); // note that SVI is more stable in the presence of multicore... #ifdef CHASTE_CVODE TS_ASSERT_DELTA(final_voltage_svi_cvode[20], svi_20, 3.0); #endif //Cvode TS_ASSERT_DELTA(final_voltage_svit[20], svi_20, 3.0); // node 130 (for h=0.02) is on the y-axis (cross-fibre direction), ICI CV is slower TS_ASSERT_DELTA(mesh.GetNode(130)->rGetLocation()[0], 0.0, 1e-9); TS_ASSERT_DELTA(mesh.GetNode(130)->rGetLocation()[1], 0.1, 1e-9); TS_ASSERT_DELTA(final_voltage_ici[130], ici_130, 1.0); TS_ASSERT_DELTA(final_voltage_svi[130], svi_130, 0.2); #ifdef CHASTE_CVODE TS_ASSERT_DELTA(final_voltage_svi_cvode[130], svi_130, 0.3); // different CVODE versions = slightly different answer! #endif //cvode TS_ASSERT_DELTA(final_voltage_svit[130], svi_130, 0.2); }
void TestConductionVelocityConvergesFasterWithSvi1d() { double h[3] = {0.001,0.01,0.02}; unsigned probe_node_index[3] = {300, 30, 15}; unsigned number_of_nodes[3] = {1001, 101, 51}; std::vector<double> conduction_vel_ici(3); std::vector<double> conduction_vel_svi(3); ReplicatableVector final_voltage_ici; ReplicatableVector final_voltage_svi; //HeartConfig::Instance()->SetUseRelativeTolerance(1e-8); HeartConfig::Instance()->SetSimulationDuration(4.0); //ms HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.01, 0.01, 0.01); for (unsigned i=0; i<3; i++) { // ICI - ionic current interpolation - the default { DistributedTetrahedralMesh<1,1> mesh; mesh.ConstructRegularSlabMesh(h[i], 1.0); TS_ASSERT_EQUALS(mesh.GetNumNodes(), number_of_nodes[i]); //Double check (for later) that the indexing is as expected if (mesh.GetDistributedVectorFactory()->IsGlobalIndexLocal( probe_node_index[i] )) { TS_ASSERT_DELTA(mesh.GetNode( probe_node_index[i] )->rGetLocation()[0], 0.3, 1e-8); } std::stringstream output_dir; output_dir << "MonodomainIci_" << h[i]; HeartConfig::Instance()->SetOutputDirectory(output_dir.str()); HeartConfig::Instance()->SetOutputFilenamePrefix("results"); // need to have this for i=1,2 cases!! HeartConfig::Instance()->SetUseStateVariableInterpolation(false); BlockCellFactory<1> cell_factory; MonodomainProblem<1> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); final_voltage_ici.ReplicatePetscVector(monodomain_problem.GetSolution()); //// see #1633 //// end time needs to be increased for these (say, to 7ms) // Hdf5DataReader simulation_data(OutputFileHandler::GetChasteTestOutputDirectory() + output_dir.str(), // "results", false); // PropagationPropertiesCalculator ppc(&simulation_data); // unsigned node_at_0_04 = (unsigned)round(0.04/h[i]); // unsigned node_at_0_40 = (unsigned)round(0.40/h[i]); // assert(fabs(mesh.GetNode(node_at_0_04)->rGetLocation()[0]-0.04)<1e-6); // assert(fabs(mesh.GetNode(node_at_0_40)->rGetLocation()[0]-0.40)<1e-6); // conduction_vel_ici[i] = ppc.CalculateConductionVelocity(node_at_0_04,node_at_0_40,0.36); // std::cout << "conduction_vel_ici = " << conduction_vel_ici[i] << "\n"; } // SVI - state variable interpolation { DistributedTetrahedralMesh<1,1> mesh; mesh.ConstructRegularSlabMesh(h[i], 1.0); //Double check (for later) that the indexing is as expected if (mesh.GetDistributedVectorFactory()->IsGlobalIndexLocal( probe_node_index[i] )) { TS_ASSERT_DELTA(mesh.GetNode( probe_node_index[i] )->rGetLocation()[0], 0.3, 1e-8); } std::stringstream output_dir; output_dir << "MonodomainSvi_" << h[i]; HeartConfig::Instance()->SetOutputDirectory(output_dir.str()); HeartConfig::Instance()->SetOutputFilenamePrefix("results"); HeartConfig::Instance()->SetUseStateVariableInterpolation(); BlockCellFactory<1> cell_factory; MonodomainProblem<1> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); final_voltage_svi.ReplicatePetscVector(monodomain_problem.GetSolution()); // Hdf5DataReader simulation_data(OutputFileHandler::GetChasteTestOutputDirectory() + output_dir.str(), // "results", false); // PropagationPropertiesCalculator ppc(&simulation_data); // unsigned node_at_0_04 = (unsigned)round(0.04/h[i]); // unsigned node_at_0_40 = (unsigned)round(0.40/h[i]); // assert(fabs(mesh.GetNode(node_at_0_04)->rGetLocation()[0]-0.04)<1e-6); // assert(fabs(mesh.GetNode(node_at_0_40)->rGetLocation()[0]-0.40)<1e-6); // conduction_vel_svi[i] = ppc.CalculateConductionVelocity(node_at_0_04,node_at_0_40,0.36); // std::cout << "conduction_vel_svi = " << conduction_vel_svi[i] << "\n"; } if (i==0) // finest mesh { for (unsigned j=0; j<final_voltage_ici.GetSize(); j++) { // visually checked they agree at this mesh resolution, and chosen tolerance from results TS_ASSERT_DELTA(final_voltage_ici[j], final_voltage_svi[j], 0.3); if (final_voltage_ici[j]>-80) { // shouldn't be exactly equal, as long as away from resting potential TS_ASSERT_DIFFERS(final_voltage_ici[j], final_voltage_svi[j]); } } double ici_voltage_at_0_03_finest_mesh = final_voltage_ici[ probe_node_index[i] ]; double svi_voltage_at_0_03_finest_mesh = final_voltage_svi[ probe_node_index[i] ]; TS_ASSERT_DELTA(svi_voltage_at_0_03_finest_mesh, 11.0067, 2e-4); //hardcoded value from fine svi TS_ASSERT_DELTA(ici_voltage_at_0_03_finest_mesh, 11.0067, 1.2e-1); //hardcoded value from fine svi } else if (i==1) { double ici_voltage_at_0_03_middle_mesh = final_voltage_ici[ probe_node_index[i] ]; double svi_voltage_at_0_03_middle_mesh = final_voltage_svi[ probe_node_index[i] ]; // ICI conduction velocity > SVI conduction velocity // and both should be greater than CV on finesh mesh TS_ASSERT_DELTA(ici_voltage_at_0_03_middle_mesh, 19.8924, 1e-3); TS_ASSERT_DELTA(svi_voltage_at_0_03_middle_mesh, 14.9579, 1e-3); } else { double ici_voltage_at_0_03_coarse_mesh = final_voltage_ici[ probe_node_index[i] ]; double svi_voltage_at_0_03_coarse_mesh = final_voltage_svi[ probe_node_index[i] ]; // ICI conduction velocity even greater than SVI conduction // velocity TS_ASSERT_DELTA(ici_voltage_at_0_03_coarse_mesh, 24.4938, 1e-3); TS_ASSERT_DELTA(svi_voltage_at_0_03_coarse_mesh, 17.3131, 1e-3); } } }
// Solve on a 2D 1mm by 1mm mesh (space step = 0.1mm), stimulating the left // edge. void TestMonodomainFitzHughNagumoWithEdgeStimulus( void ) throw (Exception) { HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(0.01, 0.01)); HeartConfig::Instance()->SetSimulationDuration(1.2); //ms HeartConfig::Instance()->SetMeshFileName("mesh/test/data/2D_0_to_1mm_400_elements"); HeartConfig::Instance()->SetOutputDirectory("FhnWithEdgeStimulus"); HeartConfig::Instance()->SetOutputFilenamePrefix("MonodomainFhn_2dWithEdgeStimulus"); FhnEdgeStimulusCellFactory cell_factory; // using the criss-cross mesh so wave propagates properly MonodomainProblem<2> monodomain_problem( &cell_factory ); monodomain_problem.Initialise(); HeartConfig::Instance()->SetSurfaceAreaToVolumeRatio(1.0); HeartConfig::Instance()->SetCapacitance(1.0); monodomain_problem.Solve(); /* * Test the top right node against the right one in the 1D case, * comparing voltage, and then test all the nodes on the right hand * side of the square against the top right one, comparing voltage. */ bool need_initialisation = true; double probe_voltage=0.0; DistributedVector voltage = monodomain_problem.GetSolutionDistributedVector(); need_initialisation = true; // Test the RHS of the mesh for (DistributedVector::Iterator node_index = voltage.Begin(); node_index != voltage.End(); ++node_index) { if (monodomain_problem.rGetMesh().GetNode(node_index.Global)->GetPoint()[0] == 0.1) { // x = 0 is where the stimulus has been applied // x = 0.1cm is the other end of the mesh and where we want to // to test the value of the nodes if (need_initialisation) { probe_voltage = voltage[node_index]; need_initialisation = false; } else { // Tests the final voltages for all the RHS edge nodes // are close to each other. // This works as we are using the 'criss-cross' mesh, // the voltages would vary more with a mesh with all the // triangles aligned in the same direction. TS_ASSERT_DELTA(voltage[node_index], probe_voltage, 2e-4); } TS_ASSERT_DELTA(voltage[node_index], 0.139426, 2e-3); } } }
void TestMonodomain3d() throw(Exception) { /* HOW_TO_TAG Cardiac/Problem definition * Generate a slab (cuboid) mesh rather than read a mesh in, and pass it to solver */ /* We will auto-generate a mesh this time, and pass it in, rather than * provide a mesh file name. This is how to generate a cuboid mesh with * a given spatial stepsize h */ TetrahedralMesh<3,3> mesh; double h=0.02; mesh.ConstructRegularSlabMesh(h, 0.8 /*length*/, 0.3 /*width*/, 0.3 /*depth*/); /* (In 2D the call is identical, but without the depth parameter). * * EMPTYLINE * * Set the simulation duration, etc, and create an instance of the cell factory. * One thing that should be noted for monodomain problems, the ''intracellular * conductivity'' is used as the monodomain effective conductivity (not a * harmonic mean of intra and extracellular conductivities). So if you want to * alter the monodomain conductivity call * `HeartConfig::Instance()->SetIntracellularConductivities` */ HeartConfig::Instance()->SetSimulationDuration(5); //ms HeartConfig::Instance()->SetOutputDirectory("Monodomain3dExample"); HeartConfig::Instance()->SetOutputFilenamePrefix("results"); HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.005, 0.01, 0.1); BenchmarkCellFactory cell_factory; /* Now we declare the problem class, `MonodomainProblem<3>` instead of `BidomainProblem<2>`. * The interface for both is the same. */ MonodomainProblem<3> monodomain_problem( &cell_factory ); /* If a mesh-file-name hasn't been set using `HeartConfig`, we have to pass in * a mesh using the `SetMesh` method (must be called before `Initialise`). */ monodomain_problem.SetMesh(&mesh); /* By default data for all nodes is output, but for big simulations, sometimes this * might not be required, and the action potential only at certain nodes required. * The following code shows how to output the results at the first, middle and last * nodes, for example. (The output is written to the HDF5 file; regular visualisation output * will be turned off. HDF5 files can be read using Matlab). We are not using this in this * simulation however (hence the boolean being set to false). */ bool partial_output = false; if(partial_output) { std::vector<unsigned> nodes_to_be_output; nodes_to_be_output.push_back(0); nodes_to_be_output.push_back((unsigned)round( (mesh.GetNumNodes()-1)/2 )); nodes_to_be_output.push_back(mesh.GetNumNodes()-1); monodomain_problem.SetOutputNodes(nodes_to_be_output); } /* `SetWriteInfo` is a useful method that means that the min/max voltage is * printed as the simulation runs (useful for verifying that cells are stimulated * and the wave propagating, for example) (although note scons does buffer output * before printing to screen) */ monodomain_problem.SetWriteInfo(); /* Finally, call `Initialise` and `Solve` as before */ monodomain_problem.Initialise(); monodomain_problem.Solve(); /* This part is just to check nothing has accidentally been changed in this example */ ReplicatableVector voltage(monodomain_problem.GetSolution()); TS_ASSERT_DELTA(voltage[0], 34.9032, 1e-2); }
void TestMonodomainTissueUsingPurkinjeCellFactory() throw(Exception) { HeartConfig::Instance()->Reset(); TrianglesMeshReader<2,2> reader("mesh/test/data/mixed_dimension_meshes/2D_0_to_1mm_200_elements"); MixedDimensionMesh<2,2> mixed_mesh; mixed_mesh.ConstructFromMeshReader(reader); PurkinjeCellFactory cell_factory; cell_factory.SetMesh(&mixed_mesh); MonodomainTissue<2> tissue( &cell_factory ); TS_ASSERT(tissue.HasPurkinje()); TS_ASSERT_EQUALS(tissue.rGetPurkinjeCellsDistributed().size(), tissue.rGetCellsDistributed().size()); TS_ASSERT_EQUALS(tissue.rGetPurkinjeIionicCacheReplicated().GetSize(), tissue.rGetIionicCacheReplicated().GetSize()); for (AbstractTetrahedralMesh<2,2>::NodeIterator current_node = mixed_mesh.GetNodeIteratorBegin(); current_node != mixed_mesh.GetNodeIteratorEnd(); ++current_node) { unsigned global_index = current_node->GetIndex(); AbstractCardiacCellInterface* p_purkinje_cell = tissue.GetPurkinjeCell(global_index); double y = current_node->rGetLocation()[1]; // cable nodes are on y=0.05 (we don't test by index because indices may be permuted in parallel). if( fabs(y-0.05) < 1e-8 ) { TS_ASSERT(dynamic_cast<CellDiFrancescoNoble1985FromCellML*>(p_purkinje_cell) != NULL); } else { TS_ASSERT(dynamic_cast<FakeBathCell*>(p_purkinje_cell) != NULL); } TS_ASSERT_EQUALS(tissue.rGetPurkinjeCellsDistributed()[global_index-mixed_mesh.GetDistributedVectorFactory()->GetLow()], p_purkinje_cell); } // Test archiving too FileFinder archive_dir("monodomain_tissue_purkinje_archive", RelativeTo::ChasteTestOutput); std::string archive_file = "monodomain_tissue_purkinje.arch"; { // Save to archive ArchiveOpener<boost::archive::text_oarchive, std::ofstream> arch_opener(archive_dir, archive_file); boost::archive::text_oarchive* p_arch = arch_opener.GetCommonArchive(); // Make sure at least one Purkinje cell has a non-initial-condition state variable to compare if (mixed_mesh.GetDistributedVectorFactory()->IsGlobalIndexLocal(0u)) { tissue.GetPurkinjeCell(0u)->SetVoltage(1234.5); } AbstractCardiacTissue<2>* const p_archive_tissue = &tissue; (*p_arch) << p_archive_tissue; } { // Load from archive and compare ArchiveOpener<boost::archive::text_iarchive, std::ifstream> arch_opener(archive_dir, archive_file); boost::archive::text_iarchive* p_arch = arch_opener.GetCommonArchive(); AbstractCardiacTissue<2>* p_tissue; (*p_arch) >> p_tissue; TS_ASSERT(p_tissue->HasPurkinje()); TS_ASSERT_EQUALS(p_tissue->rGetPurkinjeCellsDistributed().size(), tissue.rGetPurkinjeCellsDistributed().size()); TS_ASSERT_EQUALS(p_tissue->rGetPurkinjeIionicCacheReplicated().GetSize(), tissue.rGetPurkinjeIionicCacheReplicated().GetSize()); for (AbstractTetrahedralMesh<2,2>::NodeIterator current_node = p_tissue->mpMesh->GetNodeIteratorBegin(); current_node != p_tissue->mpMesh->GetNodeIteratorEnd(); ++current_node) { unsigned global_index = current_node->GetIndex(); AbstractCardiacCellInterface* p_purkinje_cell = p_tissue->GetPurkinjeCell(global_index); double y = current_node->rGetLocation()[1]; // cable nodes are on y=0.05 (we don't test by index because indices may be permuted in parallel). if( fabs(y-0.05) < 1e-8 ) { TS_ASSERT(dynamic_cast<CellDiFrancescoNoble1985FromCellML*>(p_purkinje_cell) != NULL); } else { TS_ASSERT(dynamic_cast<FakeBathCell*>(p_purkinje_cell) != NULL); } TS_ASSERT_EQUALS(p_purkinje_cell->GetVoltage(), tissue.GetPurkinjeCell(global_index)->GetVoltage()); TS_ASSERT_EQUALS(p_tissue->rGetPurkinjeCellsDistributed()[global_index-p_tissue->mpMesh->GetDistributedVectorFactory()->GetLow()], p_purkinje_cell); } delete p_tissue; } const std::string migration_archive_dir("TestMonodomainTissue/purkinje_migration_archive"); { // Save via MonodomainProblem so we can migrate // Note: from Chaste release 3.1 onward we no longer support Boost 1.33. // The earliest version of Boost supported in 1.34 // Run the test with b=_hostconfig,boost=1-34_5 to save /* scons b=_hostconfig,boost=1-34_5 ts=heart/test/monodomain/TestMonodomainTissue.hpp * */ MonodomainProblem<2> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mixed_mesh); monodomain_problem.Initialise(); TS_ASSERT(monodomain_problem.GetMonodomainTissue()->HasPurkinje()); CardiacSimulationArchiver<MonodomainProblem<2> >::Save(monodomain_problem, migration_archive_dir); } TS_ASSERT_EQUALS(tissue.rGetPurkinjeIionicCacheReplicated().GetSize(), 121u); }
// like TestOperatorSplittingMonodomainSolver but much finer mesh and smaller // dt so can check for proper convergence void TestConvergenceOnFineMesh() throw(Exception) { ReplicatableVector final_voltage_normal; ReplicatableVector final_voltage_operator_splitting; HeartConfig::Instance()->SetSimulationDuration(4.0); //ms HeartConfig::Instance()->SetOutputFilenamePrefix("results"); double h = 0.001; // very fine // Normal { HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.001, 0.001, 0.1); // very small timesteps TetrahedralMesh<1,1> mesh; mesh.ConstructRegularSlabMesh(h, 1.0); HeartConfig::Instance()->SetOutputDirectory("MonodomainCompareWithOperatorSplitting_normal"); BlockCellFactory<1> cell_factory; MonodomainProblem<1> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); final_voltage_normal.ReplicatePetscVector(monodomain_problem.GetSolution()); } // Operator splitting { HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.001, 0.002, 0.1); // very small timesteps - use pde-dt = 2 ode-dt // as effective_ode_dt = min{pde_dt/2, ode_dt} in // our operator splitting implementation TetrahedralMesh<1,1> mesh; mesh.ConstructRegularSlabMesh(h, 1.0); HeartConfig::Instance()->SetOutputDirectory("MonodomainCompareWithOperatorSplitting_splitting"); BlockCellFactory<1> cell_factory; HeartConfig::Instance()->SetUseReactionDiffusionOperatorSplitting(); MonodomainProblem<1> monodomain_problem( &cell_factory ); monodomain_problem.SetMesh(&mesh); monodomain_problem.Initialise(); monodomain_problem.Solve(); final_voltage_operator_splitting.ReplicatePetscVector(monodomain_problem.GetSolution()); } bool some_node_depolarised = false; assert(final_voltage_normal.GetSize()==final_voltage_operator_splitting.GetSize()); for(unsigned j=0; j<final_voltage_normal.GetSize(); j++) { double tol=4.7; TS_ASSERT_DELTA(final_voltage_normal[j], final_voltage_operator_splitting[j], tol); if(final_voltage_normal[j]>-80) { // shouldn't be exactly equal, as long as away from resting potential TS_ASSERT_DIFFERS(final_voltage_normal[j], final_voltage_operator_splitting[j]); } if(final_voltage_normal[j]>0.0) { some_node_depolarised = true; } } UNUSED_OPT(some_node_depolarised); assert(some_node_depolarised); }