CardiacElectroMechanicsProblem.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #include "CardiacElectroMechanicsProblem.hpp"
00030 
00031 #include "OutputFileHandler.hpp"
00032 #include "ReplicatableVector.hpp"
00033 #include "HeartConfig.hpp"
00034 #include "LogFile.hpp"
00035 #include "ChastePoint.hpp"
00036 #include "Element.hpp"
00037 #include "BoundaryConditionsContainer.hpp"
00038 #include "AbstractDynamicLinearPdeSolver.hpp"
00039 #include "TimeStepper.hpp"
00040 #include "QuadraturePointsGroup.hpp"
00041 #include "TrianglesMeshWriter.hpp"
00042 #include "Hdf5ToMeshalyzerConverter.hpp"
00043 #include "Hdf5ToCmguiConverter.hpp"
00044 #include "MeshalyzerMeshWriter.hpp"
00045 #include "PetscTools.hpp"
00046 #include "ImplicitCardiacMechanicsSolver.hpp"
00047 #include "ExplicitCardiacMechanicsSolver.hpp"
00048 #include "CmguiDeformedSolutionsWriter.hpp"
00049 #include "VoltageInterpolaterOntoMechanicsMesh.hpp"
00050 
00051 
00052 template<unsigned DIM>
00053 void CardiacElectroMechanicsProblem<DIM>::DetermineWatchedNodes()
00054 {
00055     assert(mIsWatchedLocation);
00056 
00057     // find the nearest electrics mesh node
00058     double min_dist = DBL_MAX;
00059     unsigned node_index = UNSIGNED_UNSET;
00060     for(unsigned i=0; i<mpElectricsMesh->GetNumNodes(); i++)
00061     {
00062         double dist = norm_2(mWatchedLocation - mpElectricsMesh->GetNode(i)->rGetLocation());
00063         if(dist < min_dist)
00064         {
00065             min_dist = dist;
00066             node_index = i;
00067         }
00068     }
00069 
00070     // set up watched node, if close enough
00071     assert(node_index != UNSIGNED_UNSET); // should def have found something
00072     c_vector<double,DIM> pos = mpElectricsMesh->GetNode(node_index)->rGetLocation();
00073 
00074     if(min_dist > 1e-8)
00075     {
00076         #define COVERAGE_IGNORE
00077         std::cout << "ERROR: Could not find an electrics node very close to requested watched location - "
00078                   << "min distance was " << min_dist << " for node " << node_index
00079                   << " at location " << pos << std::flush;;
00080 
00082         //EXCEPTION("Could not find an electrics node very close to requested watched location");
00083         NEVER_REACHED;
00084         #undef COVERAGE_IGNORE
00085     }
00086     else
00087     {
00088         LOG(2,"Chosen electrics node "<<node_index<<" at location " << pos << " to be watched");
00089         mWatchedElectricsNodeIndex = node_index;
00090     }
00091 
00092     // find nearest mechanics mesh
00093     min_dist = DBL_MAX;
00094     node_index = UNSIGNED_UNSET;
00095     c_vector<double,DIM> pos_at_min;
00096 
00097     for(unsigned i=0; i<mpMechanicsMesh->GetNumNodes(); i++)
00098     {
00099         c_vector<double,DIM> position = mpMechanicsMesh->GetNode(i)->rGetLocation();
00100 
00101         double dist = norm_2(position-mWatchedLocation);
00102 
00103         if(dist < min_dist)
00104         {
00105             min_dist = dist;
00106             node_index = i;
00107             pos_at_min = position;
00108         }
00109     }
00110 
00111     // set up watched node, if close enough
00112     assert(node_index != UNSIGNED_UNSET); // should def have found something
00113 
00114     if(min_dist > 1e-8)
00115     {
00116         #define COVERAGE_IGNORE
00117         std::cout << "ERROR: Could not find a mechanics node very close to requested watched location - "
00118                   << "min distance was " << min_dist << " for node " << node_index
00119                   << " at location " << pos_at_min;
00120 
00122         //EXCEPTION("Could not find a mechanics node very close to requested watched location");
00123         NEVER_REACHED;
00124         #undef COVERAGE_IGNORE
00125     }
00126     else
00127     {
00128         LOG(2,"Chosen mechanics node "<<node_index<<" at location " << pos << " to be watched");
00129         mWatchedMechanicsNodeIndex = node_index;
00130     }
00131 
00132     OutputFileHandler handler(mOutputDirectory,false);
00133     mpWatchedLocationFile = handler.OpenOutputFile("watched.txt");
00134 }
00135 
00136 template<unsigned DIM>
00137 void CardiacElectroMechanicsProblem<DIM>::WriteWatchedLocationData(double time, Vec voltage)
00138 {
00139     assert(mIsWatchedLocation);
00140 
00141     std::vector<c_vector<double,DIM> >& deformed_position = mpMechanicsSolver->rGetDeformedPosition();
00142 
00144     ReplicatableVector voltage_replicated(voltage);
00145     double V = voltage_replicated[mWatchedElectricsNodeIndex];
00146 
00149     //     // Metadata is currently being added to CellML models and then this will be avoided by asking for Calcium.
00150     //    double Ca = mpMonodomainProblem->GetMonodomainTissue()->GetCardiacCell(mWatchedElectricsNodeIndex)->GetIntracellularCalciumConcentration();
00151 
00152     *mpWatchedLocationFile << time << " ";
00153     for(unsigned i=0; i<DIM; i++)
00154     {
00155         *mpWatchedLocationFile << deformed_position[mWatchedMechanicsNodeIndex](i) << " ";
00156     }
00157     *mpWatchedLocationFile << V << "\n";
00158     mpWatchedLocationFile->flush();
00159 }
00160 
00161 template<unsigned DIM>
00162 c_matrix<double,DIM,DIM>& CardiacElectroMechanicsProblem<DIM>::rGetModifiedConductivityTensor(unsigned elementIndex, const c_matrix<double,DIM,DIM>& rOriginalConductivity)
00163 {
00164     if(mLastModifiedConductivity.first==elementIndex)
00165     {
00166         // if already computed on this electrics element, just return conductivity
00167         return mLastModifiedConductivity.second;
00168     }
00169     else
00170     {
00171         // new element, need to compute conductivity
00172 
00173         // first get the deformation gradient for this electrics element
00174         unsigned containing_mechanics_elem = mpMeshPair->rGetCoarseElementsForFineElementCentroids()[elementIndex];
00175         c_matrix<double,DIM,DIM>& r_deformation_gradient = mDeformationGradientsForEachMechanicsElement[containing_mechanics_elem];
00176 
00177         // compute sigma_def = F^{-1} sigma_undef F^{-T}
00178         c_matrix<double,DIM,DIM> inv_F = Inverse(r_deformation_gradient);
00179         mLastModifiedConductivity.second = prod(inv_F, rOriginalConductivity);
00180         mLastModifiedConductivity.second = prod(mLastModifiedConductivity.second, trans(inv_F));
00181 
00182         // save the current element and return the tensor
00183         mLastModifiedConductivity.first = elementIndex;
00184         return mLastModifiedConductivity.second;
00185     }
00186 }
00187 
00188 
00189 //
00190 //
00192 //template<unsigned DIM>
00193 //void CardiacElectroMechanicsProblem<DIM>::SetImpactRegion(std::vector<BoundaryElement<DIM-1,DIM>*>& rImpactRegion)
00194 //{
00195 //    assert(mpImpactRegion == NULL);
00196 //    mpImpactRegion = &rImpactRegion;
00197 //}
00198 //
00200 //template<unsigned DIM>
00201 //void CardiacElectroMechanicsProblem<DIM>::ApplyImpactTractions(double time)
00202 //{
00203 //    if(mpImpactRegion==NULL)
00204 //    {
00205 //        return;
00206 //    }
00207 //
00208 //    double start_time = 10;
00209 //    double end_time = 20;
00210 //    double magnitude = 1.5;
00211 //    unsigned direction = 1;
00212 //
00213 //    c_vector<double,DIM> traction = zero_vector<double>(DIM);
00214 //    if( (time>=start_time) && (time<=end_time) )
00215 //    {
00216 //        traction(direction) = magnitude;
00217 //    }
00218 //    mImpactTractions.clear();
00219 //    mImpactTractions.resize(mpImpactRegion->size(), traction);
00220 //    mpCardiacMechSolver->SetSurfaceTractionBoundaryConditions(*mpImpactRegion, mImpactTractions);
00221 //}
00222 
00223 template<unsigned DIM>
00224 CardiacElectroMechanicsProblem<DIM>::CardiacElectroMechanicsProblem(
00225             CompressibilityType compressibilityType,
00226             TetrahedralMesh<DIM,DIM>* pElectricsMesh,
00227             QuadraticMesh<DIM>* pMechanicsMesh,
00228             AbstractCardiacCellFactory<DIM>* pCellFactory,
00229             ElectroMechanicsProblemDefinition<DIM>* pProblemDefinition,
00230             std::string outputDirectory = "")
00231       : mCompressibilityType(compressibilityType),
00232         mpCardiacMechSolver(NULL),
00233         mpMechanicsSolver(NULL),
00234         mpElectricsMesh(pElectricsMesh),
00235         mpMechanicsMesh(pMechanicsMesh),
00236         mpProblemDefinition(pProblemDefinition),
00237         mpMeshPair(NULL),
00238         mNoElectricsOutput(false),
00239         mIsWatchedLocation(false),
00240         mWatchedElectricsNodeIndex(UNSIGNED_UNSET),
00241         mWatchedMechanicsNodeIndex(UNSIGNED_UNSET)
00242 {
00243     // Do some initial set up...
00244     // However, NOTE, we don't use either the passed in meshes or the problem_definition.
00245     // These pointers are allowed to be NULL, in case a child constructor wants to set
00246     // them up (eg CardiacElectroMechProbRegularGeom).
00247     // The meshes and problem_defn are used for the first time in Initialise().
00248 
00249 
00250     // Start-up mechanics event handler..
00251     MechanicsEventHandler::Reset();
00252     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL);
00253     // disable the electric event handler, because we use a problem class but
00254     // don't call Solve, so we would have to worry about starting and ending any
00255     // events in AbstractCardiacProblem::Solve() (esp. calling EndEvent(EVERYTHING))
00256     // if we didn't disable it.
00257     HeartEventHandler::Disable();
00258 
00259     assert(HeartConfig::Instance()->GetSimulationDuration()>0.0);
00260     assert(HeartConfig::Instance()->GetPdeTimeStep()>0.0);
00261 
00262     // Create the monodomain problem.
00263     // **NOTE** WE ONLY USE THIS TO: set up the cells, get an initial condition
00264     // (voltage) vector, and get an solver. We won't ever call solve on the MonodomainProblem
00265     assert(pCellFactory != NULL);
00266     mpMonodomainProblem = new MonodomainProblem<DIM>(pCellFactory);
00267 
00268     // check whether output is required
00269     mWriteOutput = (outputDirectory!="");
00270     if(mWriteOutput)
00271     {
00272         mOutputDirectory = outputDirectory;
00273         // create the directory
00274         OutputFileHandler handler(mOutputDirectory);
00275         mDeformationOutputDirectory = mOutputDirectory + "/deformation";
00276         HeartConfig::Instance()->SetOutputDirectory(mOutputDirectory + "/electrics");
00277         HeartConfig::Instance()->SetOutputFilenamePrefix("voltage");
00278     }
00279     else
00280     {
00281         mDeformationOutputDirectory = "";
00282     }
00283 
00284     mLastModifiedConductivity.first = UNSIGNED_UNSET; //important
00285 
00286 //    mpImpactRegion=NULL;
00287 }
00288 
00289 template<unsigned DIM>
00290 CardiacElectroMechanicsProblem<DIM>::~CardiacElectroMechanicsProblem()
00291 {
00292     // NOTE if SetWatchedLocation but not Initialise has been called, mpWatchedLocationFile
00293     // will be uninitialised and using it will cause a seg fault. Hence the mpMechanicsMesh!=NULL
00294     // it is true if Initialise has been called.
00295     if(mIsWatchedLocation && mpMechanicsMesh)
00296     {
00297         mpWatchedLocationFile->close();
00298     }
00299 
00300     delete mpMonodomainProblem;
00301     delete mpCardiacMechSolver;
00302     delete mpMeshPair;
00303 
00304     LogFile::Close();
00305 }
00306 
00307 template<unsigned DIM>
00308 void CardiacElectroMechanicsProblem<DIM>::Initialise()
00309 {
00310     assert(mpElectricsMesh!=NULL);
00311     assert(mpMechanicsMesh!=NULL);
00312     assert(mpProblemDefinition!=NULL);
00313 
00314     assert(mpCardiacMechSolver==NULL);
00315 
00316 
00317     mNumElecTimestepsPerMechTimestep = (unsigned) floor((mpProblemDefinition->GetMechanicsSolveTimestep()/HeartConfig::Instance()->GetPdeTimeStep())+0.5);
00318     if(fabs(mNumElecTimestepsPerMechTimestep*HeartConfig::Instance()->GetPdeTimeStep() - mpProblemDefinition->GetMechanicsSolveTimestep()) > 1e-6)
00319     {
00320         EXCEPTION("Electrics PDE timestep does not divide mechanics solve timestep");
00321     }
00322 
00323     // Create the Logfile (note we have to do this after the output dir has been
00324     // created, else the log file might get cleaned away
00325     std::string log_dir = mOutputDirectory; // just the TESTOUTPUT dir if mOutputDir="";
00326     LogFile::Instance()->Set(2, mOutputDirectory);
00327     LogFile::Instance()->WriteHeader("Electromechanics");
00328     LOG(2, DIM << "d Implicit CardiacElectroMechanics Simulation:");
00329     LOG(2, "End time = " << HeartConfig::Instance()->GetSimulationDuration() << ", electrics time step = " << HeartConfig::Instance()->GetPdeTimeStep() << ", mechanics timestep = " << mpProblemDefinition->GetMechanicsSolveTimestep() << "\n");
00330     LOG(2, "Contraction model ode timestep " << mpProblemDefinition->GetContractionModelOdeTimestep());
00331     LOG(2, "Output is written to " << mOutputDirectory << "/[deformation/electrics]");
00332 
00333     LOG(2, "Electrics mesh has " << mpElectricsMesh->GetNumNodes() << " nodes");
00334     LOG(2, "Mechanics mesh has " << mpMechanicsMesh->GetNumNodes() << " nodes");
00335 
00336     LOG(2, "Initialising..");
00337 
00338 
00339     if(mIsWatchedLocation)
00340     {
00341         DetermineWatchedNodes();
00342     }
00343 
00344     // initialise monodomain problem
00345     mpMonodomainProblem->SetMesh(mpElectricsMesh);
00346     mpMonodomainProblem->Initialise();
00347 
00348 
00349     if(mCompressibilityType==INCOMPRESSIBLE)
00350     {
00351         // NOTE: if adding to this, do below as well for compressible option.
00352         switch(mpProblemDefinition->GetContractionModel())
00353         {
00354             case NASH2004:
00355                 // Create an EXPLICIT, INCOMPRESSIBLE solver
00356                 mpCardiacMechSolver = new ExplicitCardiacMechanicsSolver<IncompressibleNonlinearElasticitySolver<DIM>,DIM>(
00357                         mpProblemDefinition->GetContractionModel(),*mpMechanicsMesh,*mpProblemDefinition,mDeformationOutputDirectory);
00358             break;
00359 
00360             case KERCHOFFS2003:
00361             case NHS:
00362                 // Create an IMPLICIT, INCOMPRESSIBLE solver
00363                 mpCardiacMechSolver = new ImplicitCardiacMechanicsSolver<IncompressibleNonlinearElasticitySolver<DIM>,DIM>(
00364                         mpProblemDefinition->GetContractionModel(),*mpMechanicsMesh,*mpProblemDefinition,mDeformationOutputDirectory);
00365             break;
00366             default:
00367                 EXCEPTION("Invalid contraction model, options are: NASH2004, KERCHOFFS2003 or NHS");
00368         }
00369     }
00370     else
00371     {
00372         // repeat above with Compressible solver rather than incompressible - not the neatest but avoids having to template
00373         // this class.
00374         switch(mpProblemDefinition->GetContractionModel())
00375         {
00376             case NASH2004:
00377                 // Create an EXPLICIT, COMPRESSIBLE solver
00378                 mpCardiacMechSolver = new ExplicitCardiacMechanicsSolver<CompressibleNonlinearElasticitySolver<DIM>,DIM>(
00379                         mpProblemDefinition->GetContractionModel(),*mpMechanicsMesh,*mpProblemDefinition,mDeformationOutputDirectory);
00380             break;
00381 
00382             case KERCHOFFS2003:
00383             case NHS:
00384                 // Create an IMPLICIT, COMPRESSIBLE solver
00385                 mpCardiacMechSolver = new ImplicitCardiacMechanicsSolver<CompressibleNonlinearElasticitySolver<DIM>,DIM>(
00386                         mpProblemDefinition->GetContractionModel(),*mpMechanicsMesh,*mpProblemDefinition,mDeformationOutputDirectory);
00387             break;
00388             default:
00389                 EXCEPTION("Invalid contraction model, options are: NASH2004, KERCHOFFS2003 or NHS");
00390         }
00391     }
00392 
00393     mpMechanicsSolver = dynamic_cast<AbstractNonlinearElasticitySolver<DIM>*>(mpCardiacMechSolver);
00394     assert(mpMechanicsSolver);
00395 
00396     if(mpProblemDefinition->ReadFibreSheetDirectionsFromFile())
00397     {
00398        mpCardiacMechSolver->SetVariableFibreSheetDirections(mpProblemDefinition->GetFibreSheetDirectionsFile(),
00399                                                             mpProblemDefinition->GetFibreSheetDirectionsDefinedPerQuadraturePoint());
00400     }
00401 
00402     // set up mesh pair and determine the fine mesh elements and corresponding weights for each
00403     // quadrature point in the coarse mesh
00404     mpMeshPair = new FineCoarseMeshPair<DIM>(*mpElectricsMesh, *mpMechanicsMesh);
00405     mpMeshPair->SetUpBoxesOnFineMesh();
00406     mpMeshPair->ComputeFineElementsAndWeightsForCoarseQuadPoints(*(mpCardiacMechSolver->GetQuadratureRule()), false);
00407     mpMeshPair->DeleteFineBoxCollection();
00408 
00409 
00410 
00411     if(mpProblemDefinition->GetDeformationAffectsConductivity() || mpProblemDefinition->GetDeformationAffectsCellModels())
00412     {
00413         mpMeshPair->SetUpBoxesOnCoarseMesh();
00414     }
00415 
00416 
00417     if(mpProblemDefinition->GetDeformationAffectsCellModels() || mpProblemDefinition->GetDeformationAffectsConductivity())
00418     {
00419         // initialise the stretches saved for each mechanics element
00420         mStretchesForEachMechanicsElement.resize(mpMechanicsMesh->GetNumElements(),1.0);
00421 
00422         // initialise the store of the F in each mechanics element (one constant value of F) in each
00423         mDeformationGradientsForEachMechanicsElement.resize(mpMechanicsMesh->GetNumElements(),identity_matrix<double>(DIM));
00424     }
00425 
00426 
00427     if(mpProblemDefinition->GetDeformationAffectsCellModels())
00428     {
00429         // compute the coarse elements which contain each fine node -- for transferring stretch from
00430         // mechanics solve electrics cell models
00431         mpMeshPair->ComputeCoarseElementsForFineNodes(false);
00432 
00433     }
00434 
00435     if(mpProblemDefinition->GetDeformationAffectsConductivity())
00436     {
00437         // compute the coarse elements which contain each fine element centroid -- for transferring F from
00438         // mechanics solve to electrics mesh elements
00439         mpMeshPair->ComputeCoarseElementsForFineElementCentroids(false);
00440 
00441         // tell the abstract tissue class that the conductivities need to be modified, passing in this class
00442         // (which is of type AbstractConductivityModifier)
00443         mpMonodomainProblem->GetMonodomainTissue()->SetConductivityModifier(this);
00444     }
00445 
00446     if(mWriteOutput)
00447     {
00448         TrianglesMeshWriter<DIM,DIM> mesh_writer(mOutputDirectory,"electrics_mesh",false);
00449         mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00450     }
00451 }
00452 
00453 template<unsigned DIM>
00454 void CardiacElectroMechanicsProblem<DIM>::Solve()
00455 {
00456     // initialise the meshes and mechanics solver
00457     if(mpCardiacMechSolver==NULL)
00458     {
00459         Initialise();
00460     }
00461 
00462     mpProblemDefinition->Validate();
00463 
00464     boost::shared_ptr<BoundaryConditionsContainer<DIM,DIM,1> > p_bcc(new BoundaryConditionsContainer<DIM,DIM,1>);
00465     p_bcc->DefineZeroNeumannOnMeshBoundary(mpElectricsMesh, 0);
00466     mpMonodomainProblem->SetBoundaryConditionsContainer(p_bcc);
00467 
00468     // get an electrics solver from the problem. Note that we don't call
00469     // Solve() on the CardiacProblem class, we do the looping here.
00470     AbstractDynamicLinearPdeSolver<DIM,DIM,1>* p_electrics_solver
00471        = mpMonodomainProblem->CreateSolver();
00472 
00473     // set up initial voltage etc
00474     Vec voltage=NULL; //This will be set and used later
00475     Vec initial_voltage = mpMonodomainProblem->CreateInitialCondition();
00476 
00477     unsigned num_quad_points = mpCardiacMechSolver->GetTotalNumQuadPoints();
00478     std::vector<double> interpolated_calcium_concs(num_quad_points, 0.0);
00479     std::vector<double> interpolated_voltages(num_quad_points, 0.0);
00480 
00481     // write the initial position
00482     unsigned counter = 0;
00483 
00484     TimeStepper stepper(0.0, HeartConfig::Instance()->GetSimulationDuration(), mpProblemDefinition->GetMechanicsSolveTimestep());
00485 
00486     CmguiDeformedSolutionsWriter<DIM>* p_cmgui_writer = NULL;
00487 
00488     unsigned mech_writer_counter = 0;
00489     if (mWriteOutput)
00490     {
00491         mpMechanicsSolver->SetWriteOutput();
00492         mpMechanicsSolver->WriteCurrentDeformation("solution",mech_writer_counter);
00493 
00494         p_cmgui_writer = new CmguiDeformedSolutionsWriter<DIM>(mOutputDirectory+"/deformation/cmgui",
00495                                                                "solution",
00496                                                                *(this->mpMechanicsMesh),
00497                                                                WRITE_QUADRATIC_MESH);
00498         std::vector<std::string> fields;
00499         fields.push_back("V");
00500         p_cmgui_writer->SetAdditionalFieldNames(fields);
00501         p_cmgui_writer->WriteInitialMesh();
00502 
00503 
00504         if(!mNoElectricsOutput)
00505         {
00506             // the writer inside monodomain problem uses the printing timestep
00507             // inside HeartConfig to estimate total number of timesteps, so make
00508             // sure this is set to what we will use.
00509             HeartConfig::Instance()->SetPrintingTimeStep(mpProblemDefinition->GetMechanicsSolveTimestep());
00510             mpMonodomainProblem->InitialiseWriter();
00511             mpMonodomainProblem->WriteOneStep(stepper.GetTime(), initial_voltage);
00512         }
00513 
00514         if(mIsWatchedLocation)
00515         {
00516             WriteWatchedLocationData(stepper.GetTime(), initial_voltage);
00517         }
00518     }
00519 
00520     PrepareForSolve();
00521 
00522     while (!stepper.IsTimeAtEnd())
00523     {
00524         LOG(2, "\nCurrent time = " << stepper.GetTime());
00525         #ifdef MECH_VERBOSE // defined in AbstractNonlinearElasticitySolver
00526         // also output time to screen as newton solve information will be output
00527         std::cout << "\n\n ** Current time = " << stepper.GetTime() << "\n";
00528         #endif
00529 
00536         if(mpProblemDefinition->GetDeformationAffectsCellModels() || mpProblemDefinition->GetDeformationAffectsConductivity())
00537         {
00538             //  Determine the stretch and deformation gradient on each element.
00539             //
00540             //  If mpProblemDefinition->GetDeformationAffectsCellModels()==true:
00541             //  Stretch will be passed to the cell models.
00542             //
00543             //  If mpProblemDefinition->GetDeformationAffectsConductivity()==true:
00544             //  The deformation gradient needs to be set up but does not need to be passed to the tissue
00545             //  so that F is used to compute the conductivity. Instead this is
00546             //  done through the line "mpMonodomainProblem->GetMonodomainTissue()->SetConductivityModifier(this);" line above, which means
00547             //  rGetModifiedConductivityTensor() will be called on this class by the tissue, which then uses the F
00548 
00549             mpCardiacMechSolver->ComputeDeformationGradientAndStretchInEachElement(mDeformationGradientsForEachMechanicsElement, mStretchesForEachMechanicsElement);
00550         }
00551 
00552         if( mpProblemDefinition->GetDeformationAffectsCellModels() )
00553         {
00554             //  Set the stretches on each of the cell models
00555             for(unsigned global_index = mpElectricsMesh->GetDistributedVectorFactory()->GetLow();
00556                          global_index < mpElectricsMesh->GetDistributedVectorFactory()->GetHigh();
00557                          global_index++)
00558             {
00559                 unsigned containing_elem = mpMeshPair->rGetCoarseElementsForFineNodes()[global_index];
00560                 double stretch = mStretchesForEachMechanicsElement[containing_elem];
00561                 mpMonodomainProblem->GetTissue()->GetCardiacCell(global_index)->SetStretch(stretch);
00562             }
00563         }
00564 
00565         p_electrics_solver->SetTimeStep(HeartConfig::Instance()->GetPdeTimeStep());
00566 
00572         LOG(2, "  Solving electrics");
00573         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::NON_MECH);
00574         for(unsigned i=0; i<mNumElecTimestepsPerMechTimestep; i++)
00575         {
00576             double current_time = stepper.GetTime() + i*HeartConfig::Instance()->GetPdeTimeStep();
00577             double next_time = stepper.GetTime() + (i+1)*HeartConfig::Instance()->GetPdeTimeStep();
00578 
00579             // solve the electrics
00580             p_electrics_solver->SetTimes(current_time, next_time);
00581             p_electrics_solver->SetInitialCondition( initial_voltage );
00582 
00583             voltage = p_electrics_solver->Solve();
00584 
00585             PetscReal min_voltage, max_voltage;
00586             VecMax(voltage,PETSC_NULL,&max_voltage); //the second param is where the index would be returned
00587             VecMin(voltage,PETSC_NULL,&min_voltage);
00588             if(i==0)
00589             {
00590                 LOG(2, "  minimum and maximum voltage is " << min_voltage <<", "<<max_voltage);
00591             }
00592             else if(i==1)
00593             {
00594                 LOG(2, "  ..");
00595             }
00596 
00597             VecDestroy(initial_voltage);
00598             initial_voltage = voltage;
00599         }
00600 
00601         if(mpProblemDefinition->GetDeformationAffectsConductivity())
00602         {
00603             p_electrics_solver->SetMatrixIsNotAssembled();
00604         }
00605 
00606 
00612 
00613         // compute Ca_I at each quad point (by interpolation, using the info on which
00614         // electrics element the quad point is in. Then set Ca_I on the mechanics solver
00615         LOG(2, "  Interpolating Ca_I and voltage");
00616 
00617         ReplicatableVector voltage_repl(voltage);
00618 
00619         // collect all the calcium concentrations (from the cells, which are
00620         // distributed) in one (replicated) vector
00621         ReplicatableVector calcium_repl(mpElectricsMesh->GetNumNodes());
00622         PetscInt lo, hi;
00623         VecGetOwnershipRange(voltage, &lo, &hi);
00624         for(int index=lo; index<hi; index++)
00625         {
00626             calcium_repl[index] = mpMonodomainProblem->GetTissue()->GetCardiacCell(index)->GetIntracellularCalciumConcentration();
00627         }
00628         PetscTools::Barrier();
00629         calcium_repl.Replicate(lo,hi);
00630 
00631 
00632         for(unsigned i=0; i<mpMeshPair->rGetElementsAndWeights().size(); i++)
00633         {
00634             double interpolated_CaI = 0;
00635             double interpolated_voltage = 0;
00636 
00637             Element<DIM,DIM>& element = *(mpElectricsMesh->GetElement(mpMeshPair->rGetElementsAndWeights()[i].ElementNum));
00638             for(unsigned node_index = 0; node_index<element.GetNumNodes(); node_index++)
00639             {
00640                 unsigned global_node_index = element.GetNodeGlobalIndex(node_index);
00641                 double CaI_at_node =  calcium_repl[global_node_index];
00642                 interpolated_CaI += CaI_at_node*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00643                 interpolated_voltage += voltage_repl[global_node_index]*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00644             }
00645 
00646             interpolated_calcium_concs[i] = interpolated_CaI;
00647             interpolated_voltages[i] = interpolated_voltage;
00648         }
00649 
00650         LOG(2, "  Setting Ca_I. max value = " << Max(interpolated_calcium_concs));
00651 
00652         // NOTE IF NHS: HERE WE SHOULD PERHAPS CHECK WHETHER THE CELL MODELS HAVE Ca_Trop
00653         // AND UPDATE FROM NHS TO CELL_MODEL, BUT NOT SURE HOW TO DO THIS.. (esp for implicit)
00654 
00655         // set [Ca], V, t
00656         mpCardiacMechSolver->SetCalciumAndVoltage(interpolated_calcium_concs, interpolated_voltages);
00657         MechanicsEventHandler::EndEvent(MechanicsEventHandler::NON_MECH);
00658 
00659 
00660 
00666         LOG(2, "  Solving mechanics ");
00667         mpMechanicsSolver->SetWriteOutput(false);
00668 
00669         // make sure the mechanics solver knows the current time (in case
00670         // the traction say is time-dependent).
00671         mpMechanicsSolver->SetCurrentTime(stepper.GetTime());
00672 
00674 //        ApplyImpactTractions(stepper.GetTime());
00675 
00676         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL_MECH);
00677         mpCardiacMechSolver->Solve(stepper.GetTime(), stepper.GetNextTime(), mpProblemDefinition->GetContractionModelOdeTimestep());
00678         MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL_MECH);
00679 
00680         LOG(2, "    Number of newton iterations = " << mpMechanicsSolver->GetNumNewtonIterations());
00681 
00682          // update the current time
00683         stepper.AdvanceOneTimeStep();
00684         counter++;
00685 
00686 
00692         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::OUTPUT);
00693         if(mWriteOutput && (counter%WRITE_EVERY_NTH_TIME==0))
00694         {
00695             LOG(2, "  Writing output");
00696             // write deformed position
00697             mech_writer_counter++;
00698             mpMechanicsSolver->SetWriteOutput();
00699             mpMechanicsSolver->WriteCurrentDeformation("solution",mech_writer_counter);
00700 
00701             p_cmgui_writer->WriteDeformationPositions(rGetDeformedPosition(), counter);
00702 
00703             if(!mNoElectricsOutput)
00704             {
00705                 mpMonodomainProblem->mpWriter->AdvanceAlongUnlimitedDimension();
00706                 mpMonodomainProblem->WriteOneStep(stepper.GetTime(), voltage);
00707             }
00708 
00709             if(mIsWatchedLocation)
00710             {
00711                 WriteWatchedLocationData(stepper.GetTime(), voltage);
00712             }
00713             OnEndOfTimeStep(counter);
00714         }
00715         MechanicsEventHandler::EndEvent(MechanicsEventHandler::OUTPUT);
00716 
00717         // write the total elapsed time..
00718         LogFile::Instance()->WriteElapsedTime("  ");
00719     }
00720 
00721 
00722 
00723     if ((mWriteOutput) && (!mNoElectricsOutput))
00724     {
00725         HeartConfig::Instance()->Reset();
00726         mpMonodomainProblem->mpWriter->Close();
00727         delete mpMonodomainProblem->mpWriter;
00728 
00729         std::string input_dir = mOutputDirectory+"/electrics";
00730 
00731         // Convert simulation data to meshalyzer format - commented
00732         std::string config_directory = HeartConfig::Instance()->GetOutputDirectory();
00733         HeartConfig::Instance()->SetOutputDirectory(input_dir);
00734 
00735         //Hdf5ToMeshalyzerConverter<DIM,DIM> meshalyzer_converter(input_dir, "voltage", mpElectricsMesh);
00736 
00737         // convert output to CMGUI format
00738         Hdf5ToCmguiConverter<DIM,DIM> cmgui_converter(input_dir,"voltage",mpElectricsMesh);
00739 
00740         // Write mesh in a suitable form for meshalyzer
00741         //std::string output_directory =  mOutputDirectory + "/electrics/output";
00742         // Write the mesh
00743         //MeshalyzerMeshWriter<DIM,DIM> mesh_writer(output_directory, "mesh", false);
00744         //mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00745         // Write the parameters out
00746         //HeartConfig::Instance()->Write();
00747 
00748         // interpolate the electrical data onto the mechanics mesh nodes and write CMGUI...
00749         // Note: this calculates the data on ALL nodes of the mechanics mesh (incl internal,
00750         // non-vertex ones), which won't be used if linear CMGUI visualisation
00751         // of the mechanics solution is used.
00752         VoltageInterpolaterOntoMechanicsMesh<DIM> cnverter(*mpElectricsMesh,*mpMechanicsMesh,input_dir,"voltage");
00753 
00754         // reset to the default value
00755         HeartConfig::Instance()->SetOutputDirectory(config_directory);
00756     }
00757 
00758     if(p_cmgui_writer)
00759     {
00760         if(mNoElectricsOutput)
00761         {
00762             p_cmgui_writer->WriteCmguiScript();
00763         }
00764         else
00765         {
00766             p_cmgui_writer->WriteCmguiScript("../../electrics/cmgui_output/voltage_mechanics_mesh");
00767         }
00768         delete p_cmgui_writer;
00769     }
00770     VecDestroy(voltage);
00771     delete p_electrics_solver;
00772 
00773     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL);
00774 }
00775 
00776 
00777 
00778 template<unsigned DIM>
00779 double CardiacElectroMechanicsProblem<DIM>::Max(std::vector<double>& vec)
00780 {
00781     double max = -1e200;
00782     for(unsigned i=0; i<vec.size(); i++)
00783     {
00784         if(vec[i]>max) max=vec[i];
00785     }
00786     return max;
00787 }
00788 
00789 template<unsigned DIM>
00790 void CardiacElectroMechanicsProblem<DIM>::SetNoElectricsOutput()
00791 {
00792     mNoElectricsOutput = true;
00793 }
00794 
00795 template<unsigned DIM>
00796 void CardiacElectroMechanicsProblem<DIM>::SetWatchedPosition(c_vector<double,DIM> watchedLocation)
00797 {
00798     mIsWatchedLocation = true;
00799     mWatchedLocation = watchedLocation;
00800 }
00801 
00802 
00803 template<unsigned DIM>
00804 std::vector<c_vector<double,DIM> >& CardiacElectroMechanicsProblem<DIM>::rGetDeformedPosition()
00805 {
00806     return mpMechanicsSolver->rGetDeformedPosition();
00807 }
00808 
00809 
00810 
00811 
00813 // Explicit instantiation
00815 
00816 //note: 1d incompressible material doesn't make sense
00817 template class CardiacElectroMechanicsProblem<2>;
00818 template class CardiacElectroMechanicsProblem<3>;
Generated on Thu Dec 22 13:00:06 2011 for Chaste by  doxygen 1.6.3