Chaste Release::3.1
CardiacElectroMechanicsProblem.cpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #include "CardiacElectroMechanicsProblem.hpp"
00037 
00038 #include "OutputFileHandler.hpp"
00039 #include "ReplicatableVector.hpp"
00040 #include "HeartConfig.hpp"
00041 #include "LogFile.hpp"
00042 #include "ChastePoint.hpp"
00043 #include "Element.hpp"
00044 #include "BoundaryConditionsContainer.hpp"
00045 #include "AbstractDynamicLinearPdeSolver.hpp"
00046 #include "TimeStepper.hpp"
00047 #include "QuadraturePointsGroup.hpp"
00048 #include "TrianglesMeshWriter.hpp"
00049 #include "Hdf5ToMeshalyzerConverter.hpp"
00050 #include "Hdf5ToCmguiConverter.hpp"
00051 #include "MeshalyzerMeshWriter.hpp"
00052 #include "PetscTools.hpp"
00053 #include "ImplicitCardiacMechanicsSolver.hpp"
00054 #include "ExplicitCardiacMechanicsSolver.hpp"
00055 #include "CmguiDeformedSolutionsWriter.hpp"
00056 #include "VoltageInterpolaterOntoMechanicsMesh.hpp"
00057 #include "BidomainProblem.hpp"
00058 
00059 
00060 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00061 void CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::DetermineWatchedNodes()
00062 {
00063     assert(mIsWatchedLocation);
00064 
00065     // find the nearest electrics mesh node
00066     double min_dist = DBL_MAX;
00067     unsigned node_index = UNSIGNED_UNSET;
00068     for(unsigned i=0; i<mpElectricsMesh->GetNumNodes(); i++)
00069     {
00070         double dist = norm_2(mWatchedLocation - mpElectricsMesh->GetNode(i)->rGetLocation());
00071         if(dist < min_dist)
00072         {
00073             min_dist = dist;
00074             node_index = i;
00075         }
00076     }
00077 
00078     // set up watched node, if close enough
00079     assert(node_index != UNSIGNED_UNSET); // should def have found something
00080     c_vector<double,DIM> pos = mpElectricsMesh->GetNode(node_index)->rGetLocation();
00081 
00082     if(min_dist > 1e-8)
00083     {
00084         #define COVERAGE_IGNORE
00085         std::cout << "ERROR: Could not find an electrics node very close to requested watched location - "
00086                   << "min distance was " << min_dist << " for node " << node_index
00087                   << " at location " << pos << std::flush;;
00088 
00090         //EXCEPTION("Could not find an electrics node very close to requested watched location");
00091         NEVER_REACHED;
00092         #undef COVERAGE_IGNORE
00093     }
00094     else
00095     {
00096         LOG(2,"Chosen electrics node "<<node_index<<" at location " << pos << " to be watched");
00097         mWatchedElectricsNodeIndex = node_index;
00098     }
00099 
00100     // find nearest mechanics mesh
00101     min_dist = DBL_MAX;
00102     node_index = UNSIGNED_UNSET;
00103     c_vector<double,DIM> pos_at_min;
00104 
00105     for(unsigned i=0; i<mpMechanicsMesh->GetNumNodes(); i++)
00106     {
00107         c_vector<double,DIM> position = mpMechanicsMesh->GetNode(i)->rGetLocation();
00108 
00109         double dist = norm_2(position-mWatchedLocation);
00110 
00111         if(dist < min_dist)
00112         {
00113             min_dist = dist;
00114             node_index = i;
00115             pos_at_min = position;
00116         }
00117     }
00118 
00119     // set up watched node, if close enough
00120     assert(node_index != UNSIGNED_UNSET); // should def have found something
00121 
00122     if(min_dist > 1e-8)
00123     {
00124         #define COVERAGE_IGNORE
00125         std::cout << "ERROR: Could not find a mechanics node very close to requested watched location - "
00126                   << "min distance was " << min_dist << " for node " << node_index
00127                   << " at location " << pos_at_min;
00128 
00130         //EXCEPTION("Could not find a mechanics node very close to requested watched location");
00131         NEVER_REACHED;
00132         #undef COVERAGE_IGNORE
00133     }
00134     else
00135     {
00136         LOG(2,"Chosen mechanics node "<<node_index<<" at location " << pos << " to be watched");
00137         mWatchedMechanicsNodeIndex = node_index;
00138     }
00139 
00140     OutputFileHandler handler(mOutputDirectory,false);
00141     mpWatchedLocationFile = handler.OpenOutputFile("watched.txt");
00142 }
00143 
00144 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00145 void CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::WriteWatchedLocationData(double time, Vec voltage)
00146 {
00147     assert(mIsWatchedLocation);
00148 
00149     std::vector<c_vector<double,DIM> >& deformed_position = mpMechanicsSolver->rGetDeformedPosition();
00150 
00152     ReplicatableVector voltage_replicated(voltage);
00153     double V = voltage_replicated[mWatchedElectricsNodeIndex];
00154 
00157     //     // Metadata is currently being added to CellML models and then this will be avoided by asking for Calcium.
00158     //    double Ca = mpElectricsProblem->GetMonodomainTissue()->GetCardiacCell(mWatchedElectricsNodeIndex)->GetIntracellularCalciumConcentration();
00159 
00160     *mpWatchedLocationFile << time << " ";
00161     for(unsigned i=0; i<DIM; i++)
00162     {
00163         *mpWatchedLocationFile << deformed_position[mWatchedMechanicsNodeIndex](i) << " ";
00164     }
00165     *mpWatchedLocationFile << V << "\n";
00166     mpWatchedLocationFile->flush();
00167 }
00168 
00169 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00170 c_matrix<double,DIM,DIM>& CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::rGetModifiedConductivityTensor(unsigned elementIndex, const c_matrix<double,DIM,DIM>& rOriginalConductivity)
00171 {
00172     if(mLastModifiedConductivity.first==elementIndex)
00173     {
00174         // if already computed on this electrics element, just return conductivity
00175         return mLastModifiedConductivity.second;
00176     }
00177     else
00178     {
00179         // new element, need to compute conductivity
00180 
00181         // first get the deformation gradient for this electrics element
00182         unsigned containing_mechanics_elem = mpMeshPair->rGetCoarseElementsForFineElementCentroids()[elementIndex];
00183         c_matrix<double,DIM,DIM>& r_deformation_gradient = mDeformationGradientsForEachMechanicsElement[containing_mechanics_elem];
00184 
00185         // compute sigma_def = F^{-1} sigma_undef F^{-T}
00186         c_matrix<double,DIM,DIM> inv_F = Inverse(r_deformation_gradient);
00187         mLastModifiedConductivity.second = prod(inv_F, rOriginalConductivity);
00188         mLastModifiedConductivity.second = prod(mLastModifiedConductivity.second, trans(inv_F));
00189 
00190         // save the current element and return the tensor
00191         mLastModifiedConductivity.first = elementIndex;
00192         return mLastModifiedConductivity.second;
00193     }
00194 }
00195 
00196 
00200 template<unsigned DIM, unsigned PROBLEM_DIM>
00201 class CreateElectricsProblem
00202 {
00203 public:
00209     static AbstractCardiacProblem<DIM, DIM, PROBLEM_DIM>* Create(ElectricsProblemType problemType,
00210                                                                  AbstractCardiacCellFactory<DIM>* pCellFactory);
00211 };
00212 
00216 template<unsigned DIM>
00217 class CreateElectricsProblem<DIM, 1u>
00218 {
00219 public:
00225     static AbstractCardiacProblem<DIM, DIM, 1u>* Create(ElectricsProblemType problemType,
00226                                                         AbstractCardiacCellFactory<DIM>* pCellFactory)
00227     {
00228         if (problemType == MONODOMAIN)
00229         {
00230             return new MonodomainProblem<DIM>(pCellFactory);
00231         }
00232         EXCEPTION("The second template parameter should be 2 when a bidomain problem is chosen");
00233     }
00234 };
00235 
00239 template<unsigned DIM>
00240 class CreateElectricsProblem<DIM, 2u>
00241 {
00242 public:
00248     static AbstractCardiacProblem<DIM, DIM, 2u>* Create(ElectricsProblemType problemType,
00249                                                         AbstractCardiacCellFactory<DIM>* pCellFactory)
00250     {
00251         if (problemType == BIDOMAIN)
00252         {
00253             return new BidomainProblem<DIM>(pCellFactory, false);//false-> no bath
00254         }
00255         if (problemType == BIDOMAIN_WITH_BATH)
00256         {
00257             return new BidomainProblem<DIM>(pCellFactory, true);// true-> bath
00258         }
00259         EXCEPTION("The second template parameter should be 1 when a monodomain problem is chosen");
00260     }
00261 };
00262 
00263 
00264 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00265 CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::CardiacElectroMechanicsProblem(
00266             CompressibilityType compressibilityType,
00267             ElectricsProblemType electricsProblemType,
00268             TetrahedralMesh<DIM,DIM>* pElectricsMesh,
00269             QuadraticMesh<DIM>* pMechanicsMesh,
00270             AbstractCardiacCellFactory<DIM>* pCellFactory,
00271             ElectroMechanicsProblemDefinition<DIM>* pProblemDefinition,
00272             std::string outputDirectory = "")
00273       : mCompressibilityType(compressibilityType),
00274         mpCardiacMechSolver(NULL),
00275         mpMechanicsSolver(NULL),
00276         mpElectricsMesh(pElectricsMesh),
00277         mpMechanicsMesh(pMechanicsMesh),
00278         mpProblemDefinition(pProblemDefinition),
00279         mHasBath(false),
00280         mpMeshPair(NULL),
00281         mNoElectricsOutput(false),
00282         mIsWatchedLocation(false),
00283         mWatchedElectricsNodeIndex(UNSIGNED_UNSET),
00284         mWatchedMechanicsNodeIndex(UNSIGNED_UNSET),
00285         mNumTimestepsToOutputDeformationGradientsAndStress(UNSIGNED_UNSET)
00286 {
00287     // Do some initial set up...
00288     // However, NOTE, we don't use either the passed in meshes or the problem_definition.
00289     // These pointers are allowed to be NULL, in case a child constructor wants to set
00290     // them up (eg CardiacElectroMechProbRegularGeom).
00291     // The meshes and problem_defn are used for the first time in Initialise().
00292 
00293 
00294     // Start-up mechanics event handler..
00295     MechanicsEventHandler::Reset();
00296     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL);
00297     // disable the electric event handler, because we use a problem class but
00298     // don't call Solve, so we would have to worry about starting and ending any
00299     // events in AbstractCardiacProblem::Solve() (esp. calling EndEvent(EVERYTHING))
00300     // if we didn't disable it.
00301     HeartEventHandler::Disable();
00302 
00303     assert(HeartConfig::Instance()->GetSimulationDuration()>0.0);
00304     assert(HeartConfig::Instance()->GetPdeTimeStep()>0.0);
00305 
00306     // Create the monodomain problem.
00307     // **NOTE** WE ONLY USE THIS TO: set up the cells, get an initial condition
00308     // (voltage) vector, and get a solver. We won't ever call Solve on the cardiac problem class
00309     assert(pCellFactory != NULL);
00310     mpElectricsProblem = CreateElectricsProblem<DIM,ELEC_PROB_DIM>::Create(electricsProblemType, pCellFactory);
00311 
00312     if (electricsProblemType == BIDOMAIN_WITH_BATH)
00313     {
00314         mHasBath = true;
00315     }
00316     // check whether output is required
00317     mWriteOutput = (outputDirectory!="");
00318     if(mWriteOutput)
00319     {
00320         mOutputDirectory = outputDirectory;
00321         // create the directory
00322         OutputFileHandler handler(mOutputDirectory);
00323         mDeformationOutputDirectory = mOutputDirectory + "/deformation";
00324         HeartConfig::Instance()->SetOutputDirectory(mOutputDirectory + "/electrics");
00325         HeartConfig::Instance()->SetOutputFilenamePrefix("voltage");
00326     }
00327     else
00328     {
00329         mDeformationOutputDirectory = "";
00330     }
00331 
00332     mLastModifiedConductivity.first = UNSIGNED_UNSET; //important
00333 
00334 //    mpImpactRegion=NULL;
00335 }
00336 
00337 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00338 CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::~CardiacElectroMechanicsProblem()
00339 {
00340     // NOTE if SetWatchedLocation but not Initialise has been called, mpWatchedLocationFile
00341     // will be uninitialised and using it will cause a seg fault. Hence the mpMechanicsMesh!=NULL
00342     // it is true if Initialise has been called.
00343     if(mIsWatchedLocation && mpMechanicsMesh)
00344     {
00345         mpWatchedLocationFile->close();
00346     }
00347 
00348     delete mpElectricsProblem;
00349     delete mpCardiacMechSolver;
00350     delete mpMeshPair;
00351 
00352     LogFile::Close();
00353 }
00354 
00355 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00356 void CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::Initialise()
00357 {
00358     assert(mpElectricsMesh!=NULL);
00359     assert(mpMechanicsMesh!=NULL);
00360     assert(mpProblemDefinition!=NULL);
00361     assert(mpCardiacMechSolver==NULL);
00362 
00365     mNumElecTimestepsPerMechTimestep = (unsigned) floor((mpProblemDefinition->GetMechanicsSolveTimestep()/HeartConfig::Instance()->GetPdeTimeStep())+0.5);
00366     if(fabs(mNumElecTimestepsPerMechTimestep*HeartConfig::Instance()->GetPdeTimeStep() - mpProblemDefinition->GetMechanicsSolveTimestep()) > 1e-6)
00367     {
00368         EXCEPTION("Electrics PDE timestep does not divide mechanics solve timestep");
00369     }
00370 
00371     // Create the Logfile (note we have to do this after the output dir has been
00372     // created, else the log file might get cleaned away
00373     std::string log_dir = mOutputDirectory; // just the TESTOUTPUT dir if mOutputDir="";
00374     LogFile::Instance()->Set(2, mOutputDirectory);
00375     LogFile::Instance()->WriteHeader("Electromechanics");
00376     LOG(2, DIM << "d Implicit CardiacElectroMechanics Simulation:");
00377     LOG(2, "End time = " << HeartConfig::Instance()->GetSimulationDuration() << ", electrics time step = " << HeartConfig::Instance()->GetPdeTimeStep() << ", mechanics timestep = " << mpProblemDefinition->GetMechanicsSolveTimestep() << "\n");
00378     LOG(2, "Contraction model ode timestep " << mpProblemDefinition->GetContractionModelOdeTimestep());
00379     LOG(2, "Output is written to " << mOutputDirectory << "/[deformation/electrics]");
00380 
00381     LOG(2, "Electrics mesh has " << mpElectricsMesh->GetNumNodes() << " nodes");
00382     LOG(2, "Mechanics mesh has " << mpMechanicsMesh->GetNumNodes() << " nodes");
00383 
00384     LOG(2, "Initialising..");
00385 
00386 
00387     if(mIsWatchedLocation)
00388     {
00389         DetermineWatchedNodes();
00390     }
00391 
00392     // initialise electrics problem
00393     mpElectricsProblem->SetMesh(mpElectricsMesh);
00394     mpElectricsProblem->Initialise();
00395 
00396 
00397     if(mCompressibilityType==INCOMPRESSIBLE)
00398     {
00399         // NOTE: if adding to this, do below as well for compressible option.
00400         switch(mpProblemDefinition->GetContractionModel())
00401         {
00402             case NASH2004:
00403                 // Create an EXPLICIT, INCOMPRESSIBLE solver
00404                 mpCardiacMechSolver = new ExplicitCardiacMechanicsSolver<IncompressibleNonlinearElasticitySolver<DIM>,DIM>(
00405                         mpProblemDefinition->GetContractionModel(),*mpMechanicsMesh,*mpProblemDefinition,mDeformationOutputDirectory);
00406             break;
00407 
00408             case KERCHOFFS2003:
00409             case NHS:
00410                 // Create an IMPLICIT, INCOMPRESSIBLE solver
00411                 mpCardiacMechSolver = new ImplicitCardiacMechanicsSolver<IncompressibleNonlinearElasticitySolver<DIM>,DIM>(
00412                         mpProblemDefinition->GetContractionModel(),*mpMechanicsMesh,*mpProblemDefinition,mDeformationOutputDirectory);
00413             break;
00414             default:
00415                 EXCEPTION("Invalid contraction model, options are: NASH2004, KERCHOFFS2003 or NHS");
00416         }
00417     }
00418     else
00419     {
00420         // repeat above with Compressible solver rather than incompressible - not the neatest but avoids having to template
00421         // this class.
00422         switch(mpProblemDefinition->GetContractionModel())
00423         {
00424             case NASH2004:
00425                 // Create an EXPLICIT, COMPRESSIBLE solver
00426                 mpCardiacMechSolver = new ExplicitCardiacMechanicsSolver<CompressibleNonlinearElasticitySolver<DIM>,DIM>(
00427                         mpProblemDefinition->GetContractionModel(),*mpMechanicsMesh,*mpProblemDefinition,mDeformationOutputDirectory);
00428             break;
00429 
00430             case KERCHOFFS2003:
00431             case NHS:
00432                 // Create an IMPLICIT, COMPRESSIBLE solver
00433                 mpCardiacMechSolver = new ImplicitCardiacMechanicsSolver<CompressibleNonlinearElasticitySolver<DIM>,DIM>(
00434                         mpProblemDefinition->GetContractionModel(),*mpMechanicsMesh,*mpProblemDefinition,mDeformationOutputDirectory);
00435             break;
00436             default:
00437                 EXCEPTION("Invalid contraction model, options are: NASH2004, KERCHOFFS2003 or NHS");
00438         }
00439     }
00440 
00441 
00442     mpMechanicsSolver = dynamic_cast<AbstractNonlinearElasticitySolver<DIM>*>(mpCardiacMechSolver);
00443     assert(mpMechanicsSolver);
00444 
00445     // set up mesh pair and determine the fine mesh elements and corresponding weights for each
00446     // quadrature point in the coarse mesh
00447     mpMeshPair = new FineCoarseMeshPair<DIM>(*mpElectricsMesh, *mpMechanicsMesh);
00448     mpMeshPair->SetUpBoxesOnFineMesh();
00449     mpMeshPair->ComputeFineElementsAndWeightsForCoarseQuadPoints(*(mpCardiacMechSolver->GetQuadratureRule()), false);
00450     mpMeshPair->DeleteFineBoxCollection();
00451 
00452     mpCardiacMechSolver->SetFineCoarseMeshPair(mpMeshPair);
00453     mpCardiacMechSolver->Initialise();
00454 
00455     unsigned num_quad_points = mpCardiacMechSolver->GetTotalNumQuadPoints();
00456     mInterpolatedCalciumConcs.assign(num_quad_points, 0.0);
00457     mInterpolatedVoltages.assign(num_quad_points, 0.0);
00458 
00459     if(mpProblemDefinition->ReadFibreSheetDirectionsFromFile())
00460     {
00461        mpCardiacMechSolver->SetVariableFibreSheetDirections(mpProblemDefinition->GetFibreSheetDirectionsFile(),
00462                                                             mpProblemDefinition->GetFibreSheetDirectionsDefinedPerQuadraturePoint());
00463     }
00464 
00465 
00466     if(mpProblemDefinition->GetDeformationAffectsConductivity() || mpProblemDefinition->GetDeformationAffectsCellModels())
00467     {
00468         mpMeshPair->SetUpBoxesOnCoarseMesh();
00469     }
00470 
00471 
00472     if(mpProblemDefinition->GetDeformationAffectsCellModels() || mpProblemDefinition->GetDeformationAffectsConductivity())
00473     {
00474         // initialise the stretches saved for each mechanics element
00475         mStretchesForEachMechanicsElement.resize(mpMechanicsMesh->GetNumElements(),1.0);
00476 
00477         // initialise the store of the F in each mechanics element (one constant value of F) in each
00478         mDeformationGradientsForEachMechanicsElement.resize(mpMechanicsMesh->GetNumElements(),identity_matrix<double>(DIM));
00479     }
00480 
00481 
00482     if(mpProblemDefinition->GetDeformationAffectsCellModels())
00483     {
00484         // compute the coarse elements which contain each fine node -- for transferring stretch from
00485         // mechanics solve electrics cell models
00486         mpMeshPair->ComputeCoarseElementsForFineNodes(false);
00487 
00488     }
00489 
00490     if(mpProblemDefinition->GetDeformationAffectsConductivity())
00491     {
00492         // compute the coarse elements which contain each fine element centroid -- for transferring F from
00493         // mechanics solve to electrics mesh elements
00494         mpMeshPair->ComputeCoarseElementsForFineElementCentroids(false);
00495 
00496         // tell the abstract tissue class that the conductivities need to be modified, passing in this class
00497         // (which is of type AbstractConductivityModifier)
00498         mpElectricsProblem->GetTissue()->SetConductivityModifier(this);
00499     }
00500 
00501     if(mWriteOutput)
00502     {
00503         TrianglesMeshWriter<DIM,DIM> mesh_writer(mOutputDirectory,"electrics_mesh",false);
00504         mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00505     }
00506 }
00507 
00508 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00509 void CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::Solve()
00510 {
00511     // initialise the meshes and mechanics solver
00512     if(mpCardiacMechSolver==NULL)
00513     {
00514         Initialise();
00515     }
00516 
00517     bool verbose_during_solve = (    mpProblemDefinition->GetVerboseDuringSolve()
00518                                   || CommandLineArguments::Instance()->OptionExists("-mech_verbose")
00519                                   || CommandLineArguments::Instance()->OptionExists("-mech_very_verbose"));
00520 
00521 
00522     mpProblemDefinition->Validate();
00523 
00524     boost::shared_ptr<BoundaryConditionsContainer<DIM,DIM,ELEC_PROB_DIM> > p_bcc(new BoundaryConditionsContainer<DIM,DIM,ELEC_PROB_DIM>);
00525     p_bcc->DefineZeroNeumannOnMeshBoundary(mpElectricsMesh, 0);
00526     mpElectricsProblem->SetBoundaryConditionsContainer(p_bcc);
00527 
00528     // get an electrics solver from the problem. Note that we don't call
00529     // Solve() on the CardiacProblem class, we do the looping here.
00530     AbstractDynamicLinearPdeSolver<DIM,DIM,ELEC_PROB_DIM>* p_electrics_solver
00531        = mpElectricsProblem->CreateSolver();
00532 
00533     // set up initial voltage etc
00534     Vec electrics_solution=NULL; //This will be set and used later
00535     Vec calcium_data= mpElectricsMesh->GetDistributedVectorFactory()->CreateVec();
00536     Vec initial_voltage = mpElectricsProblem->CreateInitialCondition();
00537 
00538     // write the initial position
00539     unsigned counter = 0;
00540 
00541     TimeStepper stepper(0.0, HeartConfig::Instance()->GetSimulationDuration(), mpProblemDefinition->GetMechanicsSolveTimestep());
00542 
00543     CmguiDeformedSolutionsWriter<DIM>* p_cmgui_writer = NULL;
00544 
00545     std::vector<std::string> variable_names;
00546 
00547     if (mWriteOutput)
00548     {
00549         mpMechanicsSolver->SetWriteOutput();
00550         mpMechanicsSolver->WriteCurrentSpatialSolution("undeformed","nodes");
00551 
00552         p_cmgui_writer = new CmguiDeformedSolutionsWriter<DIM>(mOutputDirectory+"/deformation/cmgui",
00553                                                                "solution",
00554                                                                *(this->mpMechanicsMesh),
00555                                                                WRITE_QUADRATIC_MESH);
00556         variable_names.push_back("V");
00557         if(ELEC_PROB_DIM==2)
00558         {
00559             variable_names.push_back("Phi_e");
00560             if (mHasBath==true)
00561             {
00562                 std::vector<std::string> regions;
00563                 regions.push_back("tissue");
00564                 regions.push_back("bath");
00565                 p_cmgui_writer->SetRegionNames(regions);
00566             }
00567         }
00568         p_cmgui_writer->SetAdditionalFieldNames(variable_names);
00569         p_cmgui_writer->WriteInitialMesh("undeformed");
00570     }
00571 
00578 
00579     LOG(2, "\nSolving for initial deformation");
00580     #define COVERAGE_IGNORE
00581     if(verbose_during_solve)
00582     {
00583         std::cout << "\n\n ** Solving for initial deformation\n";
00584     }
00585     #undef COVERAGE_IGNORE
00586 
00587     mpMechanicsSolver->SetWriteOutput(false);
00588 
00589     mpMechanicsSolver->SetCurrentTime(0.0);
00590     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL_MECH);
00591 
00592     mpMechanicsSolver->SetIncludeActiveTension(false);
00593     if(mNumTimestepsToOutputDeformationGradientsAndStress!=UNSIGNED_UNSET)
00594     {
00595         mpMechanicsSolver->SetComputeAverageStressPerElementDuringSolve(true);
00596     }
00597 
00598     unsigned total_newton_iters = 0;
00599     for(unsigned index=1; index<=mpProblemDefinition->GetNumIncrementsForInitialDeformation(); index++)
00600     {
00601         #define COVERAGE_IGNORE
00602         if(verbose_during_solve)
00603         {
00604             std::cout << "    Increment " << index << " of " << mpProblemDefinition->GetNumIncrementsForInitialDeformation() << "\n";
00605         }
00606         #undef COVERAGE_IGNORE
00607 
00608         if(mpProblemDefinition->GetTractionBoundaryConditionType()==PRESSURE_ON_DEFORMED)
00609         {
00610             mpProblemDefinition->SetPressureScaling(((double)index)/mpProblemDefinition->GetNumIncrementsForInitialDeformation());
00611         }
00612         mpMechanicsSolver->Solve();
00613 
00614         total_newton_iters += mpMechanicsSolver->GetNumNewtonIterations();
00615     }
00616 
00617     mpMechanicsSolver->SetIncludeActiveTension(true);
00618     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL_MECH);
00619     LOG(2, "    Number of newton iterations = " << total_newton_iters);
00620 
00621 
00622     unsigned mech_writer_counter = 0;
00623 
00624     if (mWriteOutput)
00625     {
00626         LOG(2, "  Writing output");
00627         mpMechanicsSolver->SetWriteOutput();
00628         mpMechanicsSolver->WriteCurrentSpatialSolution("solution","nodes",mech_writer_counter);
00629         p_cmgui_writer->WriteDeformationPositions(rGetDeformedPosition(), mech_writer_counter);
00630 
00631         if(!mNoElectricsOutput)
00632         {
00633             // the writer inside monodomain problem uses the printing timestep
00634             // inside HeartConfig to estimate total number of timesteps, so make
00635             // sure this is set to what we will use.
00636             HeartConfig::Instance()->SetPrintingTimeStep(mpProblemDefinition->GetMechanicsSolveTimestep());
00637 
00638             // When we consider archiving these simulations we will need to get a bool back from the following
00639             // command to decide whether or not the file is being extended, and hence whether to write the
00640             // initial conditions into the .h5 file.
00641             mpElectricsProblem->InitialiseWriter();
00642             mpElectricsProblem->WriteOneStep(stepper.GetTime(), initial_voltage);
00643         }
00644 
00645         if(mIsWatchedLocation)
00646         {
00647             WriteWatchedLocationData(stepper.GetTime(), initial_voltage);
00648         }
00649 
00650         if(mNumTimestepsToOutputDeformationGradientsAndStress!=UNSIGNED_UNSET)
00651         {
00652             mpMechanicsSolver->WriteCurrentDeformationGradients("deformation_gradient",mech_writer_counter);
00653             mpMechanicsSolver->WriteCurrentAverageElementStresses("second_PK",mech_writer_counter);
00654         }
00655     }
00656 
00657 
00658     PrepareForSolve();
00659 
00663 //    std::vector<double> current_solution_previous_time_step = mpMechanicsSolver->rGetCurrentSolution();
00664 //    std::vector<double> current_solution_second_last_time_step = mpMechanicsSolver->rGetCurrentSolution();
00665 //    bool first_step = true;
00666 
00667 
00668     // reset this to false, may be reset again below
00669     mpMechanicsSolver->SetComputeAverageStressPerElementDuringSolve(false);
00670 
00671     while (!stepper.IsTimeAtEnd())
00672     {
00673         LOG(2, "\nCurrent time = " << stepper.GetTime());
00674         #define COVERAGE_IGNORE
00675         if(verbose_during_solve)
00676         {
00677             // also output time to screen as newton solve information will be output
00678             std::cout << "\n\n ** Current time = " << stepper.GetTime() << "\n";
00679         }
00680         #undef COVERAGE_IGNORE
00681 
00688         if(mpProblemDefinition->GetDeformationAffectsCellModels() || mpProblemDefinition->GetDeformationAffectsConductivity())
00689         {
00690             //  Determine the stretch and deformation gradient on each element.
00691             //
00692             //  If mpProblemDefinition->GetDeformationAffectsCellModels()==true:
00693             //  Stretch will be passed to the cell models.
00694             //
00695             //  If mpProblemDefinition->GetDeformationAffectsConductivity()==true:
00696             //  The deformation gradient needs to be set up but does not need to be passed to the tissue
00697             //  so that F is used to compute the conductivity. Instead this is
00698             //  done through the line "mpElectricsProblem->GetMonodomainTissue()->SetConductivityModifier(this);" line above, which means
00699             //  rGetModifiedConductivityTensor() will be called on this class by the tissue, which then uses the F
00700 
00701             mpCardiacMechSolver->ComputeDeformationGradientAndStretchInEachElement(mDeformationGradientsForEachMechanicsElement, mStretchesForEachMechanicsElement);
00702         }
00703 
00704         if( mpProblemDefinition->GetDeformationAffectsCellModels() )
00705         {
00706             //  Set the stretches on each of the cell models
00707             for(unsigned global_index = mpElectricsMesh->GetDistributedVectorFactory()->GetLow();
00708                          global_index < mpElectricsMesh->GetDistributedVectorFactory()->GetHigh();
00709                          global_index++)
00710             {
00711                 unsigned containing_elem = mpMeshPair->rGetCoarseElementsForFineNodes()[global_index];
00712                 double stretch = mStretchesForEachMechanicsElement[containing_elem];
00713                 mpElectricsProblem->GetTissue()->GetCardiacCell(global_index)->SetStretch(stretch);
00714             }
00715         }
00716 
00717         p_electrics_solver->SetTimeStep(HeartConfig::Instance()->GetPdeTimeStep());
00718 
00724         LOG(2, "  Solving electrics");
00725         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::NON_MECH);
00726         for(unsigned i=0; i<mNumElecTimestepsPerMechTimestep; i++)
00727         {
00728             double current_time = stepper.GetTime() + i*HeartConfig::Instance()->GetPdeTimeStep();
00729             double next_time = stepper.GetTime() + (i+1)*HeartConfig::Instance()->GetPdeTimeStep();
00730 
00731             // solve the electrics
00732             p_electrics_solver->SetTimes(current_time, next_time);
00733             p_electrics_solver->SetInitialCondition( initial_voltage );
00734 
00735             electrics_solution = p_electrics_solver->Solve();
00736 
00737             PetscReal min_voltage, max_voltage;
00738             VecMax(electrics_solution,PETSC_NULL,&max_voltage); //the second param is where the index would be returned
00739             VecMin(electrics_solution,PETSC_NULL,&min_voltage);
00740             if(i==0)
00741             {
00742                 LOG(2, "  minimum and maximum voltage is " << min_voltage <<", "<<max_voltage);
00743             }
00744             else if(i==1)
00745             {
00746                 LOG(2, "  ..");
00747             }
00748 
00749             PetscTools::Destroy(initial_voltage);
00750             initial_voltage = electrics_solution;
00751         }
00752 
00753         if(mpProblemDefinition->GetDeformationAffectsConductivity())
00754         {
00755             p_electrics_solver->SetMatrixIsNotAssembled();
00756         }
00757 
00758 
00764 
00765         // compute Ca_I at each quad point (by interpolation, using the info on which
00766         // electrics element the quad point is in. Then set Ca_I on the mechanics solver
00767         LOG(2, "  Interpolating Ca_I and voltage");
00768 
00769         //Collect the distributed calcium data into one Vec to be later replicated
00770         for(unsigned node_index = 0; node_index<mpElectricsMesh->GetNumNodes(); node_index++)
00771         {
00772             if (mpElectricsMesh->GetDistributedVectorFactory()->IsGlobalIndexLocal(node_index))
00773             {
00774                 double calcium_value = mpElectricsProblem->GetTissue()->GetCardiacCell(node_index)->GetIntracellularCalciumConcentration();
00775                 VecSetValue(calcium_data, node_index ,calcium_value, INSERT_VALUES);
00776             }
00777         }
00778         PetscTools::Barrier();//not sure this is needed
00779 
00780         //Replicate electrics solution and calcium (replication is inside this constructor of ReplicatableVector)
00781         ReplicatableVector electrics_solution_repl(electrics_solution);//size=(number of electrics nodes)*ELEC_PROB_DIM
00782         ReplicatableVector calcium_repl(calcium_data);//size = number of electrics nodes
00783 
00784         //interpolate values onto mechanics mesh
00785         for(unsigned i=0; i<mpMeshPair->rGetElementsAndWeights().size(); i++)
00786         {
00787             double interpolated_CaI = 0;
00788             double interpolated_voltage = 0;
00789 
00790             Element<DIM,DIM>& element = *(mpElectricsMesh->GetElement(mpMeshPair->rGetElementsAndWeights()[i].ElementNum));
00791 
00792             for(unsigned node_index = 0; node_index<element.GetNumNodes(); node_index++)
00793             {
00794                 unsigned global_index = element.GetNodeGlobalIndex(node_index);
00795                 double CaI_at_node =  calcium_repl[global_index];
00796                 interpolated_CaI += CaI_at_node*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00797                 //the following line assumes interleaved solution for ELEC_PROB_DIM>1 (e.g, [Vm_0, phi_e_0, Vm1, phi_e_1...])
00798                 interpolated_voltage += electrics_solution_repl[global_index*ELEC_PROB_DIM]*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00799             }
00800 
00801             assert(i<mInterpolatedCalciumConcs.size());
00802             assert(i<mInterpolatedVoltages.size());
00803             mInterpolatedCalciumConcs[i] = interpolated_CaI;
00804             mInterpolatedVoltages[i] = interpolated_voltage;
00805         }
00806 
00807         LOG(2, "  Setting Ca_I. max value = " << Max(mInterpolatedCalciumConcs));
00808 
00809         // NOTE IF NHS: HERE WE SHOULD PERHAPS CHECK WHETHER THE CELL MODELS HAVE Ca_Trop
00810         // AND UPDATE FROM NHS TO CELL_MODEL, BUT NOT SURE HOW TO DO THIS.. (esp for implicit)
00811 
00812         // set [Ca], V, t
00813         mpCardiacMechSolver->SetCalciumAndVoltage(mInterpolatedCalciumConcs, mInterpolatedVoltages);
00814         MechanicsEventHandler::EndEvent(MechanicsEventHandler::NON_MECH);
00815 
00816 
00822         LOG(2, "  Solving mechanics ");
00823         mpMechanicsSolver->SetWriteOutput(false);
00824 
00825         // make sure the mechanics solver knows the current time (in case
00826         // the traction say is time-dependent).
00827         mpMechanicsSolver->SetCurrentTime(stepper.GetTime());
00828 
00829         // see if we will need to output stresses at the end of this timestep
00830         if(    mNumTimestepsToOutputDeformationGradientsAndStress!=UNSIGNED_UNSET
00831             && (counter+1)%mNumTimestepsToOutputDeformationGradientsAndStress == 0 )
00832         {
00833             mpMechanicsSolver->SetComputeAverageStressPerElementDuringSolve(true);
00834         }
00835 
00838 //        for(unsigned i=0; i<mpMechanicsSolver->rGetCurrentSolution().size(); i++)
00839 //        {
00840 //            double current = mpMechanicsSolver->rGetCurrentSolution()[i];
00841 //            double previous = current_solution_previous_time_step[i];
00842 //            double second_last = current_solution_second_last_time_step[i];
00843 //            //double guess = 2*current - previous;
00844 //            double guess = 3*current - 3*previous + second_last;
00845 //
00846 //            if(!first_step)
00847 //            {
00848 //                current_solution_second_last_time_step[i] = current_solution_previous_time_step[i];
00849 //            }
00850 //            current_solution_previous_time_step[i] = mpMechanicsSolver->rGetCurrentSolution()[i];
00851 //            mpMechanicsSolver->rGetCurrentSolution()[i] = guess;
00852 //        }
00853 //        first_step = false;
00854 
00855         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL_MECH);
00856         mpCardiacMechSolver->Solve(stepper.GetTime(), stepper.GetNextTime(), mpProblemDefinition->GetContractionModelOdeTimestep());
00857         MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL_MECH);
00858 
00859         LOG(2, "    Number of newton iterations = " << mpMechanicsSolver->GetNumNewtonIterations());
00860 
00861          // update the current time
00862         stepper.AdvanceOneTimeStep();
00863         counter++;
00864 
00865 
00871         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::OUTPUT);
00872         if(mWriteOutput && (counter%WRITE_EVERY_NTH_TIME==0))
00873         {
00874             LOG(2, "  Writing output");
00875             // write deformed position
00876             mech_writer_counter++;
00877             mpMechanicsSolver->SetWriteOutput();
00878             mpMechanicsSolver->WriteCurrentSpatialSolution("solution","nodes",mech_writer_counter);
00879 
00880             p_cmgui_writer->WriteDeformationPositions(rGetDeformedPosition(), counter);
00881 
00882             if(!mNoElectricsOutput)
00883             {
00884                 mpElectricsProblem->mpWriter->AdvanceAlongUnlimitedDimension();
00885                 mpElectricsProblem->WriteOneStep(stepper.GetTime(), electrics_solution);
00886             }
00887 
00888             if(mIsWatchedLocation)
00889             {
00890                 WriteWatchedLocationData(stepper.GetTime(), electrics_solution);
00891             }
00892             OnEndOfTimeStep(counter);
00893 
00894             if(mNumTimestepsToOutputDeformationGradientsAndStress!=UNSIGNED_UNSET && counter%mNumTimestepsToOutputDeformationGradientsAndStress==0)
00895             {
00896                 mpMechanicsSolver->WriteCurrentDeformationGradients("deformation_gradient",mech_writer_counter);
00897                 mpMechanicsSolver->WriteCurrentAverageElementStresses("second_PK",mech_writer_counter);
00898             }
00899             mpMechanicsSolver->SetComputeAverageStressPerElementDuringSolve(false);
00900         }
00901         MechanicsEventHandler::EndEvent(MechanicsEventHandler::OUTPUT);
00902 
00903         // write the total elapsed time..
00904         LogFile::Instance()->WriteElapsedTime("  ");
00905     }
00906 
00907 
00908 
00909     if ((mWriteOutput) && (!mNoElectricsOutput))
00910     {
00911         HeartConfig::Instance()->Reset();
00912         mpElectricsProblem->mpWriter->Close();
00913         delete mpElectricsProblem->mpWriter;
00914 
00915         std::string input_dir = mOutputDirectory+"/electrics";
00916 
00917         // Convert simulation data to meshalyzer format - commented
00918         std::string config_directory = HeartConfig::Instance()->GetOutputDirectory();
00919         HeartConfig::Instance()->SetOutputDirectory(input_dir);
00920 
00921         //Hdf5ToMeshalyzerConverter<DIM,DIM> meshalyzer_converter(input_dir, "voltage", mpElectricsMesh);
00922 
00923         // convert output to CMGUI format
00924         Hdf5ToCmguiConverter<DIM,DIM> cmgui_converter(input_dir,"voltage",mpElectricsMesh, mHasBath);
00925 
00926         // Write mesh in a suitable form for meshalyzer
00927         //std::string output_directory =  mOutputDirectory + "/electrics/output";
00928         // Write the mesh
00929         //MeshalyzerMeshWriter<DIM,DIM> mesh_writer(output_directory, "mesh", false);
00930         //mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00931         // Write the parameters out
00932         //HeartConfig::Instance()->Write();
00933 
00934         // interpolate the electrical data onto the mechanics mesh nodes and write CMGUI...
00935         // Note: this calculates the data on ALL nodes of the mechanics mesh (incl internal,
00936         // non-vertex ones), which won't be used if linear CMGUI visualisation
00937         // of the mechanics solution is used.
00938         VoltageInterpolaterOntoMechanicsMesh<DIM> converter(*mpElectricsMesh,*mpMechanicsMesh,variable_names,input_dir,"voltage");
00939 
00940         // reset to the default value
00941         HeartConfig::Instance()->SetOutputDirectory(config_directory);
00942     }
00943 
00944     if(p_cmgui_writer)
00945     {
00946         if(mNoElectricsOutput)
00947         {
00948             p_cmgui_writer->WriteCmguiScript("","undeformed");
00949         }
00950         else
00951         {
00952             p_cmgui_writer->WriteCmguiScript("../../electrics/cmgui_output/voltage_mechanics_mesh","undeformed");
00953         }
00954         delete p_cmgui_writer;
00955     }
00956     PetscTools::Destroy(electrics_solution);
00957     PetscTools::Destroy(calcium_data);
00958     delete p_electrics_solver;
00959 
00960     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL);
00961 }
00962 
00963 
00964 
00965 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00966 double CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::Max(std::vector<double>& vec)
00967 {
00968     double max = -1e200;
00969     for(unsigned i=0; i<vec.size(); i++)
00970     {
00971         if(vec[i]>max) max=vec[i];
00972     }
00973     return max;
00974 }
00975 
00976 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00977 void CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::SetNoElectricsOutput()
00978 {
00979     mNoElectricsOutput = true;
00980 }
00981 
00982 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00983 void CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::SetWatchedPosition(c_vector<double,DIM> watchedLocation)
00984 {
00985     mIsWatchedLocation = true;
00986     mWatchedLocation = watchedLocation;
00987 }
00988 
00989 template<unsigned DIM, unsigned ELEC_PROB_DIM>
00990 void CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::SetOutputDeformationGradientsAndStress(double timeStep)
00991 {
00992     mNumTimestepsToOutputDeformationGradientsAndStress = (unsigned) floor((timeStep/mpProblemDefinition->GetMechanicsSolveTimestep())+0.5);
00993     if(fabs(mNumTimestepsToOutputDeformationGradientsAndStress*mpProblemDefinition->GetMechanicsSolveTimestep() - timeStep) > 1e-6)
00994     {
00995         EXCEPTION("Timestep provided for SetOutputDeformationGradientsAndStress() is not a multiple of mechanics solve timestep");
00996     }
00997 }
00998 
00999 
01000 template<unsigned DIM, unsigned ELEC_PROB_DIM>
01001 std::vector<c_vector<double,DIM> >& CardiacElectroMechanicsProblem<DIM,ELEC_PROB_DIM>::rGetDeformedPosition()
01002 {
01003     return mpMechanicsSolver->rGetDeformedPosition();
01004 }
01005 
01006 
01007 
01008 
01010 // Explicit instantiation
01012 
01013 //note: 1d incompressible material doesn't make sense
01014 template class CardiacElectroMechanicsProblem<2,1>;
01015 template class CardiacElectroMechanicsProblem<3,1>;
01016 template class CardiacElectroMechanicsProblem<2,2>;
01017 template class CardiacElectroMechanicsProblem<3,2>;