CardiacElectroMechanicsProblem.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
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 "AbstractDynamicAssemblerMixin.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 
00047 #include "ImplicitCardiacMechanicsAssembler.hpp"
00048 #include "ExplicitCardiacMechanicsAssembler.hpp"
00049 // if including Cinv in monobidomain equations
00050 //#include "NodewiseData.hpp"
00051 
00052 #include "MooneyRivlinMaterialLaw.hpp"
00053 
00054 #include "CmguiDeformedSolutionsWriter.hpp"
00055 
00056 
00057 template<unsigned DIM>
00058 void CardiacElectroMechanicsProblem<DIM>::DetermineWatchedNodes()
00059 {
00060     assert(mIsWatchedLocation);
00061 
00062     // find the nearest electrics mesh node
00063     double min_dist = DBL_MAX;
00064     unsigned node_index = UNSIGNED_UNSET;
00065     for(unsigned i=0; i<mpElectricsMesh->GetNumNodes(); i++)
00066     {
00067         double dist = norm_2(mWatchedLocation - mpElectricsMesh->GetNode(i)->rGetLocation());
00068         if(dist < min_dist)
00069         {
00070             min_dist = dist;
00071             node_index = i;
00072         }
00073     }
00074 
00075     // set up watched node, if close enough
00076     assert(node_index != UNSIGNED_UNSET); // should def have found something
00077     c_vector<double,DIM> pos = mpElectricsMesh->GetNode(node_index)->rGetLocation();
00078 
00079     if(min_dist > 1e-8)
00080     {
00081         #define COVERAGE_IGNORE
00082         std::cout << "ERROR: Could not find an electrics node very close to requested watched location - "
00083                   << "min distance was " << min_dist << " for node " << node_index
00084                   << " at location " << pos << std::flush;;
00085 
00087         //EXCEPTION("Could not find an electrics node very close to requested watched location");
00088         NEVER_REACHED;
00089         #undef COVERAGE_IGNORE
00090     }
00091     else
00092     {
00093         LOG(2,"Chosen electrics node "<<node_index<<" at location " << pos << " to be watched");
00094         mWatchedElectricsNodeIndex = node_index;
00095     }
00096 
00097     // find nearest mechanics mesh
00098     min_dist = DBL_MAX;
00099     node_index = UNSIGNED_UNSET;
00100     c_vector<double,DIM> pos_at_min;
00101 
00102     for(unsigned i=0; i<mpMechanicsMesh->GetNumNodes(); i++)
00103     {
00104         c_vector<double,DIM> position = mpMechanicsMesh->GetNode(i)->rGetLocation();
00105 
00106         double dist = norm_2(position-mWatchedLocation);
00107 
00108         if(dist < min_dist)
00109         {
00110             min_dist = dist;
00111             node_index = i;
00112             pos_at_min = position;
00113         }
00114     }
00115 
00116     // set up watched node, if close enough
00117     assert(node_index != UNSIGNED_UNSET); // should def have found something
00118 
00119     if(min_dist > 1e-8)
00120     {
00121         #define COVERAGE_IGNORE
00122         std::cout << "ERROR: Could not find a mechanics node very close to requested watched location - "
00123                   << "min distance was " << min_dist << " for node " << node_index
00124                   << " at location " << pos_at_min;
00125 
00127         //EXCEPTION("Could not find a mechanics node very close to requested watched location");
00128         assert(0);
00129         #undef COVERAGE_IGNORE
00130     }
00131     else
00132     {
00133         LOG(2,"Chosen mechanics node "<<node_index<<" at location " << pos << " to be watched");
00134         mWatchedMechanicsNodeIndex = node_index;
00135     }
00136 
00137     OutputFileHandler handler(mOutputDirectory);
00138     mpWatchedLocationFile = handler.OpenOutputFile("watched.txt");
00139 }
00140 
00141 template<unsigned DIM>
00142 void CardiacElectroMechanicsProblem<DIM>::WriteWatchedLocationData(double time, Vec voltage)
00143 {
00144     assert(mIsWatchedLocation);
00145 
00146     std::vector<c_vector<double,DIM> >& deformed_position = mpCardiacMechAssembler->rGetDeformedPosition();
00147 
00149     ReplicatableVector voltage_replicated(voltage);
00150     double V = voltage_replicated[mWatchedElectricsNodeIndex];
00151 
00154     //     // \todo: NOTE!!! HARDCODED state variable index - assumes Lr91.
00155     //     // Metadata is currently being added to CellML models and then this will be avoided by asking for Calcium.
00156     //    double Ca = mpMonodomainProblem->GetMonodomainPde()->GetCardiacCell(mWatchedElectricsNodeIndex)->rGetStateVariables()[3];
00157 
00158     *mpWatchedLocationFile << time << " ";
00159     for(unsigned i=0; i<DIM; i++)
00160     {
00161         *mpWatchedLocationFile << deformed_position[mWatchedMechanicsNodeIndex](i) << " ";
00162     }
00163     *mpWatchedLocationFile << V << "\n";
00164     mpWatchedLocationFile->flush();
00165 }
00166 
00167 
00168 template<unsigned DIM>
00169 CardiacElectroMechanicsProblem<DIM>::CardiacElectroMechanicsProblem(
00170             ContractionModel contractionModel,
00171             TetrahedralMesh<DIM,DIM>* pElectricsMesh,
00172             QuadraticMesh<DIM>* pMechanicsMesh,
00173             std::vector<unsigned> fixedMechanicsNodes,
00174             AbstractCardiacCellFactory<DIM>* pCellFactory,
00175             double endTime,
00176             double electricsPdeTimeStep,
00177             unsigned numElecTimeStepsPerMechTimestep,
00178             double contractionModelOdeTimeStep,
00179             std::string outputDirectory = "") :
00180         mpMeshPair(NULL)
00181 {
00182     // Start-up mechanics event handler..
00183     MechanicsEventHandler::Reset();
00184     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL);
00185     // disable the electric event handler, because we use a problem class but
00186     // don't call Solve, so we would have to worry about starting and ending any
00187     // events in AbstractCardiacProblem::Solve() (esp. calling EndEvent(EVERYTHING))
00188     // if we didn't disable it.
00189     HeartEventHandler::Disable();
00190 
00191     mContractionModel = contractionModel;
00192 
00193     // create the monodomain problem. Note the we use this to set up the cells,
00194     // get an initial condition (voltage) vector, and get an assembler. We won't
00195     // ever call solve on the MonodomainProblem
00196     assert(pCellFactory != NULL);
00197     mpMonodomainProblem = new MonodomainProblem<DIM>(pCellFactory);
00198 
00199     // save time infomation
00200     assert(endTime > 0);
00201     mEndTime = endTime;
00202 
00203     assert(electricsPdeTimeStep>0);
00204     mElectricsTimeStep = electricsPdeTimeStep;
00205 
00206     assert(numElecTimeStepsPerMechTimestep>0);
00207 
00208     mNumElecTimestepsPerMechTimestep = numElecTimeStepsPerMechTimestep;
00209 
00210     mMechanicsTimeStep = mElectricsTimeStep*mNumElecTimestepsPerMechTimestep;
00211     assert(contractionModelOdeTimeStep <= mMechanicsTimeStep+1e-14);
00212     mContractionModelOdeTimeStep = contractionModelOdeTimeStep;
00213 
00214     // check whether output is required
00215     mWriteOutput = (outputDirectory!="");
00216     if(mWriteOutput)
00217     {
00218         mOutputDirectory = outputDirectory;
00219         mDeformationOutputDirectory = mOutputDirectory + "/deformation";
00220         HeartConfig::Instance()->SetOutputDirectory(mOutputDirectory + "/electrics");
00221         HeartConfig::Instance()->SetOutputFilenamePrefix("voltage");
00222     }
00223     else
00224     {
00225         mDeformationOutputDirectory = "";
00226     }
00227     mNoElectricsOutput = false;
00228 
00229     // initialise all the pointers
00230     mpElectricsMesh = pElectricsMesh; // note these are allowed to be null, in case a child constructor wants to create them
00231     mpMechanicsMesh = pMechanicsMesh;
00232     mFixedNodes = fixedMechanicsNodes;
00233 
00234     mpCardiacMechAssembler = NULL;
00235 
00236     // Create the Logfile (note we have to do this after the output dir has been
00237     // created, else the log file might get cleaned away
00238     std::string log_dir = mOutputDirectory; // just the TESTOUTPUT dir if mOutputDir="";
00239     LogFile::Instance()->Set(2, mOutputDirectory);
00240     LogFile::Instance()->WriteHeader("Electromechanics");
00241     LOG(2, DIM << "d Implicit CardiacElectroMechanics Simulation:");
00242     LOG(2, "End time = " << mEndTime << ", electrics time step = " << mElectricsTimeStep << ", mechanics timestep = " << mMechanicsTimeStep << "\n");
00243     LOG(2, "Contraction model ode timestep " << mContractionModelOdeTimeStep);
00244     LOG(2, "Output is written to " << mOutputDirectory << "/[deformation/electrics]");
00245 #define COVERAGE_IGNORE
00247     if(mpElectricsMesh != NULL)
00248     {
00249         LOG(2, "Electrics mesh has " << mpElectricsMesh->GetNumNodes() << " nodes");
00250     }
00251     if(mpMechanicsMesh != NULL)
00252     {
00253         LOG(2, "Mechanics mesh has " << mpMechanicsMesh->GetNumNodes() << " nodes");
00254     }
00255 #undef COVERAGE_IGNORE
00256     mIsWatchedLocation = false;
00257     mWatchedElectricsNodeIndex = UNSIGNED_UNSET;
00258     mWatchedMechanicsNodeIndex = UNSIGNED_UNSET;
00259 
00260     mFibreSheetDirectionsFile = "";
00261 }
00262 
00263 template<unsigned DIM>
00264 CardiacElectroMechanicsProblem<DIM>::~CardiacElectroMechanicsProblem()
00265 {   
00271     if(mIsWatchedLocation && mpMechanicsMesh)
00272     {
00273         mpWatchedLocationFile->close();
00274     }
00275 
00276     delete mpMonodomainProblem;
00277 
00278     delete mpCardiacMechAssembler;
00279 
00280     delete mpMeshPair;
00281 
00282     LogFile::Close();
00283 }
00284 
00285 template<unsigned DIM>
00286 void CardiacElectroMechanicsProblem<DIM>::Initialise()
00287 {
00288     LOG(2, "Initialising meshes and cardiac mechanics assembler..");
00289 
00290     assert(mpElectricsMesh!=NULL);
00291     assert(mpMechanicsMesh!=NULL);
00292     assert(mpCardiacMechAssembler==NULL);
00293 
00294     if(mIsWatchedLocation)
00295     {
00296         DetermineWatchedNodes();
00297     }
00298 
00299 
00300     // initialise monodomain problem
00301     mpMonodomainProblem->SetMesh(mpElectricsMesh);
00302     HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(1.75,1.75,1.75));
00303     mpMonodomainProblem->Initialise();
00304 
00305     // Construct mechanics assembler
00306     // Here we pick the best solver for each particular contraction model. Commented out versions are
00307     // for experimentation.
00308     switch(mContractionModel)
00309     {
00310         case NASH2004:
00311             // stretch and stretch-rate independent, so should use explicit
00312             mpCardiacMechAssembler = new ExplicitCardiacMechanicsAssembler<DIM>(mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes);
00313             break;
00314         case KERCHOFFS2003:
00315             // stretch independent, so should use implicit (explicit may be unstable)
00316             mpCardiacMechAssembler = new ImplicitCardiacMechanicsAssembler<DIM>(mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes);
00317             break;
00318         case NHS:
00319             // stretch and stretch-rate independent, so should definitely use implicit
00320             mpCardiacMechAssembler = new ImplicitCardiacMechanicsAssembler<DIM>(mContractionModel,mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes);
00321             break;
00322         default:
00323             EXCEPTION("Invalid contraction model, options are: KERCHOFFS2003 or NHS");
00324     }
00325 
00326     if(mFibreSheetDirectionsFile!="")
00327     {
00328        mpCardiacMechAssembler->SetVariableFibreSheetDirections(mFibreSheetDirectionsFile);
00329     }
00330 
00331     // set up mesh pair and determine the fine mesh elements and corresponding weights for each
00332     // quadrature point in the coarse mesh
00333     mpMeshPair = new FineCoarseMeshPair<DIM>(*mpElectricsMesh, *mpMechanicsMesh);
00334     mpMeshPair->SetUpBoxesOnFineMesh();
00335     mpMeshPair->ComputeFineElementsAndWeightsForCoarseQuadPoints(*(mpCardiacMechAssembler->GetQuadratureRule()), false);
00336     mpMeshPair->DeleteBoxCollection();
00337 
00338     if(mWriteOutput)
00339     {
00340         TrianglesMeshWriter<DIM,DIM> mesh_writer(mOutputDirectory,"electrics_mesh",false);
00341         mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00342     }
00343 
00345 //    // get the assembler to compute which electrics nodes are in each mechanics mesh
00346 //    mpCardiacMechAssembler->ComputeElementsContainingNodes(mpElectricsMesh);
00347 //    assert(DIM==2);
00348 //
00349 //    NodewiseData<DIM>::Instance()->AllocateMemory(mpElectricsMesh->GetNumNodes(), 3);
00350 //    std::vector<std::vector<double> >& r_c_inverse = NodewiseData<DIM>::Instance()->rGetData();
00351 //    for(unsigned i=0; i<r_c_inverse.size(); i++)
00352 //    {
00353 //        r_c_inverse[i][0] = 1.0;
00354 //        r_c_inverse[i][1] = 0.0;
00355 //        r_c_inverse[i][2] = 1.0;
00356 //    }
00357 }
00358 
00359 template<unsigned DIM>
00360 void CardiacElectroMechanicsProblem<DIM>::Solve()
00361 {
00362 
00363     // initialise the meshes and mechanics assembler
00364     if(mpCardiacMechAssembler==NULL)
00365     {
00366         Initialise();
00367     }
00368 
00369     boost::shared_ptr<BoundaryConditionsContainer<DIM,DIM,1> > p_bcc(new BoundaryConditionsContainer<DIM,DIM,1>);
00370     p_bcc->DefineZeroNeumannOnMeshBoundary(mpElectricsMesh, 0);
00371     mpMonodomainProblem->SetBoundaryConditionsContainer(p_bcc);
00372 
00373     // get an electrics assembler from the problem. Note that we don't call
00374     // Solve() on the CardiacProblem class, we do the looping here.
00375     AbstractDynamicAssemblerMixin<DIM,DIM,1>* p_electrics_assembler
00376        = mpMonodomainProblem->CreateAssembler();
00377 
00378     // set up initial voltage etc
00379     Vec voltage=NULL; //This will be set and used later
00380     Vec initial_voltage = mpMonodomainProblem->CreateInitialCondition();
00381 
00382     unsigned num_quad_points = mpCardiacMechAssembler->GetTotalNumQuadPoints();
00383     std::vector<double> interpolated_calcium_concs(num_quad_points, 0.0);
00384     std::vector<double> interpolated_voltages(num_quad_points, 0.0);
00385 
00386     // write the initial position
00387     unsigned counter = 0;
00388 
00389     TimeStepper stepper(0.0, mEndTime, mMechanicsTimeStep);
00390 
00391     CmguiDeformedSolutionsWriter<DIM>* p_cmgui_writer = NULL;
00392 
00393     unsigned mech_writer_counter = 0;
00394     if (mWriteOutput)
00395     {
00396         mpCardiacMechAssembler->SetWriteOutput();
00397         mpCardiacMechAssembler->WriteOutput(mech_writer_counter);
00398 
00399         p_cmgui_writer = new CmguiDeformedSolutionsWriter<DIM>(mOutputDirectory+"/deformation/cmgui", "solution", *(this->mpMechanicsMesh));
00400         std::vector<std::string> fields;
00401         fields.push_back("V");
00402         p_cmgui_writer->SetAdditionalFieldNames(fields);
00403         p_cmgui_writer->WriteInitialMesh();
00404 
00405 
00406         if(!mNoElectricsOutput)
00407         {
00408             mpMonodomainProblem->InitialiseWriter();
00409             mpMonodomainProblem->WriteOneStep(stepper.GetTime(), initial_voltage);
00410         }
00411 
00412         if(mIsWatchedLocation)
00413         {
00414             WriteWatchedLocationData(stepper.GetTime(), initial_voltage);
00415         }
00416     }
00417 
00418     while (!stepper.IsTimeAtEnd())
00419     {
00420         LOG(2, "\nCurrent time = " << stepper.GetTime());
00421         #ifdef MECH_VERBOSE // defined in AbstractNonlinearElasticityAssembler
00422         // also output time to screen as newton solve information will be output
00423         std::cout << "\n\n ** Current time = " << stepper.GetTime() << "\n";
00424         #endif
00425 
00426         LOG(2, "  Solving electrics");
00427         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::NON_MECH);
00428         for(unsigned i=0; i<mNumElecTimestepsPerMechTimestep; i++)
00429         {
00430             double current_time = stepper.GetTime() + i*mElectricsTimeStep;
00431             double next_time = stepper.GetTime() + (i+1)*mElectricsTimeStep;
00432 
00433             // solve the electrics
00434             p_electrics_assembler->SetTimes(current_time, next_time, mElectricsTimeStep);
00435             p_electrics_assembler->SetInitialCondition( initial_voltage );
00436 
00437             voltage = p_electrics_assembler->Solve();
00438 
00439             PetscReal min_voltage, max_voltage;
00440             VecMax(voltage,PETSC_NULL,&max_voltage); //the second param is where the index would be returned
00441             VecMin(voltage,PETSC_NULL,&min_voltage);
00442             if(i==0)
00443             {
00444                 LOG(2, "  minimum and maximum voltage is " << min_voltage <<", "<<max_voltage);
00445             }
00446             else if(i==1)
00447             {
00448                 LOG(2, "  ..");
00449             }
00450 
00451             VecDestroy(initial_voltage);
00452             initial_voltage = voltage;
00453         }
00454 
00456 //        p_electrics_assembler->SetMatrixIsNotAssembled();
00457 
00458         // compute Ca_I at each quad point (by interpolation, using the info on which
00459         // electrics element the quad point is in. Then set Ca_I on the mechanics solver
00460         LOG(2, "  Interpolating Ca_I and voltage");
00461 
00462         ReplicatableVector voltage_repl(voltage);
00463 
00464         // collect all the calcium concentrations (from the cells, which are
00465         // distributed) in one (replicated) vector
00466         ReplicatableVector calcium_repl(mpElectricsMesh->GetNumNodes());
00467         PetscInt lo, hi;
00468         VecGetOwnershipRange(voltage, &lo, &hi);
00469         for(int index=lo; index<hi; index++)
00470         {
00471             calcium_repl[index] = mpMonodomainProblem->GetPde()->GetCardiacCell(index)->GetIntracellularCalciumConcentration();
00472         }
00473         PetscTools::Barrier();
00474         calcium_repl.Replicate(lo,hi);
00475 
00476 
00477         for(unsigned i=0; i<mpMeshPair->rGetElementsAndWeights().size(); i++)
00478         {
00479             double interpolated_CaI = 0;
00480             double interpolated_voltage = 0;
00481 
00482             Element<DIM,DIM>& element = *(mpElectricsMesh->GetElement(mpMeshPair->rGetElementsAndWeights()[i].ElementNum));
00483             for(unsigned node_index = 0; node_index<element.GetNumNodes(); node_index++)
00484             {
00485                 unsigned global_node_index = element.GetNodeGlobalIndex(node_index);
00486                 double CaI_at_node =  calcium_repl[global_node_index];
00487                 interpolated_CaI += CaI_at_node*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00488                 interpolated_voltage += voltage_repl[global_node_index]*mpMeshPair->rGetElementsAndWeights()[i].Weights(node_index);
00489             }
00490 
00491             interpolated_calcium_concs[i] = interpolated_CaI;
00492             interpolated_voltages[i] = interpolated_voltage;
00493         }
00494 
00495         LOG(2, "  Setting Ca_I. max value = " << Max(interpolated_calcium_concs));
00496 
00497         // NOTE IF NHS: HERE WE SHOULD PERHAPS CHECK WHETHER THE CELL MODELS HAVE Ca_Trop
00498         // AND UPDATE FROM NHS TO CELL_MODEL, BUT NOT SURE HOW TO DO THIS.. (esp for implicit)
00499 
00500         // set [Ca], V, t
00501         mpCardiacMechAssembler->SetCalciumAndVoltage(interpolated_calcium_concs, interpolated_voltages);
00502         MechanicsEventHandler::EndEvent(MechanicsEventHandler::NON_MECH);
00503 
00504         // solve the mechanics
00505         LOG(2, "  Solving mechanics ");
00506         //double timestep = std::min(0.01, stepper.GetNextTime()-stepper.GetTime());
00507         mpCardiacMechAssembler->SetWriteOutput(false);
00508 
00509         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL_MECH);
00510         mpCardiacMechAssembler->Solve(stepper.GetTime(), stepper.GetNextTime(), mContractionModelOdeTimeStep);
00511         MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL_MECH);
00512 
00513         LOG(2, "    Number of newton iterations = " << mpCardiacMechAssembler->GetNumNewtonIterations());
00514 
00515         // update the current time
00516         stepper.AdvanceOneTimeStep();
00517         counter++;
00518 
00519         // output the results
00520         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::OUTPUT);
00521         if(mWriteOutput && (counter%WRITE_EVERY_NTH_TIME==0))
00522         {
00523             LOG(2, "  Writing output");
00524             // write deformed position
00525             mech_writer_counter++;
00526             mpCardiacMechAssembler->SetWriteOutput();
00527             mpCardiacMechAssembler->WriteOutput(mech_writer_counter);
00528 
00529             p_cmgui_writer->WriteDeformationPositions(rGetDeformedPosition(), counter);
00530 
00531             if(!mNoElectricsOutput)
00532             {
00533                 mpMonodomainProblem->mpWriter->AdvanceAlongUnlimitedDimension();
00534                 mpMonodomainProblem->WriteOneStep(stepper.GetTime(), voltage);
00535             }
00536 
00537             if(mIsWatchedLocation)
00538             {
00539                 WriteWatchedLocationData(stepper.GetTime(), voltage);
00540             }
00541         }
00542         MechanicsEventHandler::EndEvent(MechanicsEventHandler::OUTPUT);
00543 
00545 //        // setup the Cinverse data;
00546 //        std::vector<std::vector<double> >& r_c_inverse = NodewiseData<DIM>::Instance()->rGetData();
00547 //        mpCardiacMechAssembler->CalculateCinverseAtNodes(mpElectricsMesh, r_c_inverse);
00548 //
00549 //        // write lambda
00550 //        std::stringstream file_name;
00551 //        file_name << "lambda_" << mech_writer_counter << ".dat";
00552 //        mpCardiacMechAssembler->WriteLambda(mOutputDirectory,file_name.str());
00553 
00554 
00555         // write the total elapsed time..
00556         LogFile::Instance()->WriteElapsedTime("  ");
00557     }
00558 
00559     if ((mWriteOutput) && (!mNoElectricsOutput))
00560     {
00561         HeartConfig::Instance()->Reset();
00562         mpMonodomainProblem->mpWriter->Close();
00563         delete mpMonodomainProblem->mpWriter;
00564 
00565 
00566         // Convert simulation data to Meshalyzer format
00567         //
00568         std::string input_dir = mOutputDirectory+"/electrics";
00569         std::string config_directory = HeartConfig::Instance()->GetOutputDirectory();
00570         HeartConfig::Instance()->SetOutputDirectory(input_dir);
00571         Hdf5ToMeshalyzerConverter<DIM,DIM> meshalyzer_converter(input_dir, "voltage", mpElectricsMesh);
00572         Hdf5ToCmguiConverter<DIM,DIM> cmgui_converter(input_dir,"voltage",mpElectricsMesh);
00573 
00574         // Write mesh in a suitable form for meshalyzer
00575         std::string output_directory =  mOutputDirectory + "/electrics/output";
00576         // Write the mesh
00577         MeshalyzerMeshWriter<DIM,DIM> mesh_writer(output_directory, "mesh", false);
00578         mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00579         // Write the parameters out
00580         HeartConfig::Instance()->Write();
00581 
00582         // reset to the default value
00583         HeartConfig::Instance()->SetOutputDirectory(config_directory);
00584 
00585     }
00586 
00587     if(p_cmgui_writer)
00588     {
00589         p_cmgui_writer->WriteCmguiScript();
00590         delete p_cmgui_writer;
00591     }
00592     VecDestroy(voltage);
00593     delete p_electrics_assembler;
00594 
00595     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL);
00596 }
00597 
00598 
00599 
00600 template<unsigned DIM>
00601 double CardiacElectroMechanicsProblem<DIM>::Max(std::vector<double>& vec)
00602 {
00603     double max = -1e200;
00604     for(unsigned i=0; i<vec.size(); i++)
00605     {
00606         if(vec[i]>max) max=vec[i];
00607     }
00608     return max;
00609 }
00610 
00611 template<unsigned DIM>
00612 void CardiacElectroMechanicsProblem<DIM>::SetNoElectricsOutput()
00613 {
00614     mNoElectricsOutput = true;
00615 }
00616 
00617 template<unsigned DIM>
00618 void CardiacElectroMechanicsProblem<DIM>::SetWatchedPosition(c_vector<double,DIM> watchedLocation)
00619 {
00620     mIsWatchedLocation = true;
00621     mWatchedLocation = watchedLocation;
00622 }
00623 
00624 template<unsigned DIM>
00625 void CardiacElectroMechanicsProblem<DIM>::SetVariableFibreSheetDirectionsFile(std::string fibreSheetDirectionsFile)
00626 {
00627     mFibreSheetDirectionsFile = fibreSheetDirectionsFile;
00628 }
00629 
00630 template<unsigned DIM>
00631 std::vector<c_vector<double,DIM> >& CardiacElectroMechanicsProblem<DIM>::rGetDeformedPosition()
00632 {
00633     return mpCardiacMechAssembler->rGetDeformedPosition();
00634 }
00635 
00636 
00638 // Explicit instantiation
00640 
00641 //note: 1d incompressible material doesn't make sense
00642 template class CardiacElectroMechanicsProblem<2>;
00643 template class CardiacElectroMechanicsProblem<3>;

Generated by  doxygen 1.6.2