CardiacElectroMechanicsProblem.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
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 "MeshalyzerMeshWriter.hpp"
00044 #include "PetscTools.hpp"
00045 
00046 // if including Cinv in monobidomain equations
00047 //#include "NodewiseData.hpp"
00048 
00049 
00050 template<unsigned DIM>
00051 void CardiacElectroMechanicsProblem<DIM>::DetermineWatchedNodes()
00052 {
00053     assert(mIsWatchedLocation);
00054 
00055     // find the nearest electrics mesh node
00056     double min_dist = DBL_MAX;
00057     unsigned node_index = UNSIGNED_UNSET;
00058     for(unsigned i=0; i<mpElectricsMesh->GetNumNodes(); i++)
00059     {
00060         double dist = norm_2(mWatchedLocation - mpElectricsMesh->GetNode(i)->rGetLocation());
00061         if(dist < min_dist)
00062         {
00063             min_dist = dist;
00064             node_index = i;
00065         }
00066     }
00067 
00068     // set up watched node, if close enough
00069     assert(node_index != UNSIGNED_UNSET); // should def have found something
00070     c_vector<double,DIM> pos = mpElectricsMesh->GetNode(node_index)->rGetLocation();
00071 
00072     if(min_dist > 1e-8)
00073     {
00074         #define COVERAGE_IGNORE
00075         std::cout << "ERROR: Could not find an electrics node very close to requested watched location - "
00076                   << "min distance was " << min_dist << " for node " << node_index
00077                   << " at location " << pos << std::flush;;
00078 
00080         //EXCEPTION("Could not find an electrics node very close to requested watched location");
00081         NEVER_REACHED;
00082         #undef COVERAGE_IGNORE
00083     }
00084     else
00085     {
00086         LOG_AND_COUT(2,"Chosen electrics node "<<node_index<<" at location " << pos << " to be watched");
00087         mWatchedElectricsNodeIndex = node_index;
00088     }
00089 
00090     // find nearest mechanics mesh
00091     min_dist = DBL_MAX;
00092     node_index = UNSIGNED_UNSET;
00093     c_vector<double,DIM> pos_at_min;
00094 
00095     for(unsigned i=0; i<mpMechanicsMesh->GetNumNodes(); i++)
00096     {
00097         c_vector<double,DIM> position = mpMechanicsMesh->GetNode(i)->rGetLocation();
00098 
00099         double dist = norm_2(position-mWatchedLocation);
00100 
00101         if(dist < min_dist)
00102         {
00103             min_dist = dist;
00104             node_index = i;
00105             pos_at_min = position;
00106         }
00107     }
00108 
00109     // set up watched node, if close enough
00110     assert(node_index != UNSIGNED_UNSET); // should def have found something
00111 
00112     if(min_dist > 1e-8)
00113     {
00114         #define COVERAGE_IGNORE
00115         std::cout << "ERROR: Could not find a mechanics node very close to requested watched location - "
00116                   << "min distance was " << min_dist << " for node " << node_index
00117                   << " at location " << pos_at_min;
00118 
00120         //EXCEPTION("Could not find a mechanics node very close to requested watched location");
00121         assert(0);
00122         #undef COVERAGE_IGNORE
00123     }
00124     else
00125     {
00126         LOG_AND_COUT(2,"Chosen electrics node "<<node_index<<" at location " << pos << " to be watched");
00127         mWatchedMechanicsNodeIndex = node_index;
00128     }
00129 
00130     OutputFileHandler handler(mOutputDirectory);
00131     mpWatchedLocationFile = handler.OpenOutputFile("watched.txt");
00132 }
00133 
00134 template<unsigned DIM>
00135 void CardiacElectroMechanicsProblem<DIM>::WriteWatchedLocationData(double time, Vec voltage)
00136 {
00137     assert(mIsWatchedLocation);
00138 
00139     std::vector<c_vector<double,DIM> >& deformed_position = mpCardiacMechAssembler->rGetDeformedPosition();
00140 
00142     ReplicatableVector voltage_replicated(voltage);
00143     double V=voltage_replicated[mWatchedElectricsNodeIndex];
00144 
00149     double Ca = mpMonodomainProblem->GetMonodomainPde()->GetCardiacCell(mWatchedElectricsNodeIndex)->rGetStateVariables()[3];
00150 
00151     *mpWatchedLocationFile << time << " ";
00152     for(unsigned i=0; i<DIM; i++)
00153     {
00154         *mpWatchedLocationFile << deformed_position[mWatchedMechanicsNodeIndex](i) << " ";
00155     }
00156     *mpWatchedLocationFile << V <<  " " << Ca << "\n";
00157     mpWatchedLocationFile->flush();
00158 }
00159 
00160 
00161 template<unsigned DIM>
00162 CardiacElectroMechanicsProblem<DIM>::CardiacElectroMechanicsProblem(
00163             TetrahedralMesh<DIM,DIM>* pElectricsMesh,
00164             QuadraticMesh<DIM>* pMechanicsMesh,
00165             std::vector<unsigned> fixedMechanicsNodes,
00166             AbstractCardiacCellFactory<DIM>* pCellFactory,
00167             double endTime,
00168             unsigned numElecTimeStepsPerMechTimestep,
00169             double nhsOdeTimeStep,
00170             std::string outputDirectory = "")
00171 {
00172     // Start-up mechanics event handler..
00173     MechanicsEventHandler::Reset();
00174     MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL);
00175     // disable the electric event handler, because we use a problem class but
00176     // don't call Solve, so we would have to worry about starting and ending any
00177     // events in AbstractCardiacProblem::Solve() (esp. calling EndEvent(EVERYTHING))
00178     // if we didn't disable it.
00179     HeartEventHandler::Disable();
00180 
00181     // create the monodomain problem. Note the we use this to set up the cells,
00182     // get an initial condition (voltage) vector, and get an assembler. We won't
00183     // ever call solve on the MonodomainProblem
00184     assert(pCellFactory != NULL);
00185     mpMonodomainProblem = new MonodomainProblem<DIM>(pCellFactory);
00186 
00187     // save time infomation
00188     assert(endTime > 0);
00189     mEndTime = endTime;
00190     mElectricsTimeStep = 0.01;
00191     assert(numElecTimeStepsPerMechTimestep>0);
00192 
00193     mNumElecTimestepsPerMechTimestep = numElecTimeStepsPerMechTimestep;
00194 
00195     mMechanicsTimeStep = mElectricsTimeStep*mNumElecTimestepsPerMechTimestep;
00196     assert(nhsOdeTimeStep <= mMechanicsTimeStep+1e-14);
00197     mNhsOdeTimeStep = nhsOdeTimeStep;
00198 
00199     // check whether output is required
00200     mWriteOutput = (outputDirectory!="");
00201     if(mWriteOutput)
00202     {
00203         mOutputDirectory = outputDirectory;
00204         mDeformationOutputDirectory = mOutputDirectory + "/deformation";
00205         HeartConfig::Instance()->SetOutputDirectory(mOutputDirectory + "/electrics");
00206         HeartConfig::Instance()->SetOutputFilenamePrefix("voltage");
00207     }
00208     else
00209     {
00210         mDeformationOutputDirectory = "";
00211     }
00212     mNoElectricsOutput = false;
00213 
00214     // initialise all the pointers
00215     mpElectricsMesh = pElectricsMesh; // note these are allowed to be null, in case a child constructor wants to create them
00216     mpMechanicsMesh = pMechanicsMesh;
00217     mFixedNodes = fixedMechanicsNodes;
00218 
00219     mpCardiacMechAssembler = NULL;
00220 
00221     // Create the Logfile (note we have to do this after the output dir has been
00222     // created, else the log file might get cleaned away
00223     std::string log_dir = mOutputDirectory; // just the TESTOUTPUT dir if mOutputDir="";
00224     LogFile::Instance()->Set(2, mOutputDirectory);
00225     LogFile::Instance()->WriteHeader("Electromechanics");
00226     LOG(2, DIM << "d Implicit CardiacElectroMechanics Simulation:");
00227     LOG(2, "End time = " << mEndTime << ", electrics time step = " << mElectricsTimeStep << ", mechanics timestep = " << mMechanicsTimeStep << "\n");
00228     LOG(2, "Nhs ode timestep " << mNhsOdeTimeStep);
00229     LOG(2, "Output is written to " << mOutputDirectory << "/[deformation/electrics]");
00230 #define COVERAGE_IGNORE
00232     if(mpElectricsMesh != NULL)
00233     {
00234         LOG(2, "Electrics mesh has " << mpElectricsMesh->GetNumNodes() << " nodes");
00235     }
00236     if(mpMechanicsMesh != NULL)
00237     {
00238         LOG(2, "Mechanics mesh has " << mpMechanicsMesh->GetNumNodes() << " nodes");
00239     }
00240 #undef COVERAGE_IGNORE
00241     mIsWatchedLocation = false;
00242     mWatchedElectricsNodeIndex = UNSIGNED_UNSET;
00243     mWatchedMechanicsNodeIndex = UNSIGNED_UNSET;
00244 }
00245 
00246 template<unsigned DIM>
00247 CardiacElectroMechanicsProblem<DIM>::~CardiacElectroMechanicsProblem()
00248 {   
00254     if(mIsWatchedLocation && mpMechanicsMesh)
00255     {
00256         mpWatchedLocationFile->close();
00257     }
00258 
00259     delete mpMonodomainProblem;
00260     delete mpCardiacMechAssembler;
00261 
00262     LogFile::Close();
00263 }
00264 
00265 template<unsigned DIM>
00266 void CardiacElectroMechanicsProblem<DIM>::Initialise()
00267 {
00268     LOG(2, "Initialising meshes and cardiac mechanics assembler..");
00269 
00270     assert(mpElectricsMesh!=NULL);
00271     assert(mpMechanicsMesh!=NULL);
00272     assert(mpCardiacMechAssembler==NULL);
00273 
00274     if(mIsWatchedLocation)
00275     {
00276         DetermineWatchedNodes();
00277     }
00278 
00279     // initialise monodomain problem
00280     mpMonodomainProblem->SetMesh(mpElectricsMesh);
00281 
00282     HeartConfig::Instance()->SetIntracellularConductivities(Create_c_vector(1.75,1.75,1.75));
00283 
00284     mpMonodomainProblem->Initialise();
00285 
00286     // construct mechanics assembler
00287     mpCardiacMechAssembler = new ImplicitCardiacMechanicsAssembler<DIM>(mpMechanicsMesh,mDeformationOutputDirectory,mFixedNodes);
00288 
00289     // find the element nums and weights for each gauss point in the mechanics mesh
00290     mElementAndWeightsForQuadPoints.resize(mpCardiacMechAssembler->GetTotalNumQuadPoints());
00291 
00292     // get the quad point positions in the mechanics assembler
00293     QuadraturePointsGroup<DIM> quad_point_posns(*mpMechanicsMesh, *(mpCardiacMechAssembler->GetQuadratureRule()));
00294 
00295     // find the electrics element and weight for each quad point in the mechanics mesh,
00296     // and store
00297     for(unsigned i=0; i<quad_point_posns.Size(); i++)
00298     {
00299         ChastePoint<DIM> point;
00300 
00301         for(unsigned j=0;j<DIM;j++)
00302         {
00303             point.rGetLocation()[j]=quad_point_posns.Get(i)[j];
00304         }
00305 
00306         unsigned elem_index = mpElectricsMesh->GetContainingElementIndex(point);
00307         c_vector<double,DIM+1> weight = mpElectricsMesh->GetElement(elem_index)->CalculateInterpolationWeights(point);
00308 
00309         mElementAndWeightsForQuadPoints[i].ElementNum = elem_index;
00310         mElementAndWeightsForQuadPoints[i].Weights = weight;
00311     }
00312 
00313     if(mWriteOutput)
00314     {
00315         TrianglesMeshWriter<DIM,DIM> mesh_writer(mOutputDirectory,"electrics_mesh",false);
00316         mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00317     }
00318 
00320 //    // get the assembler to compute which electrics nodes are in each mechanics mesh
00321 //    mpCardiacMechAssembler->ComputeElementsContainingNodes(mpElectricsMesh);
00322 //    assert(DIM==2);
00323 //
00324 //    NodewiseData<DIM>::Instance()->AllocateMemory(mpElectricsMesh->GetNumNodes(), 3);
00325 //    std::vector<std::vector<double> >& r_c_inverse = NodewiseData<DIM>::Instance()->rGetData();
00326 //    for(unsigned i=0; i<r_c_inverse.size(); i++)
00327 //    {
00328 //        r_c_inverse[i][0] = 1.0;
00329 //        r_c_inverse[i][1] = 0.0;
00330 //        r_c_inverse[i][2] = 1.0;
00331 //    }
00332 }
00333 
00334 template<unsigned DIM>
00335 void CardiacElectroMechanicsProblem<DIM>::Solve()
00336 {
00337     // initialise the meshes and mechanics assembler
00338     if(mpCardiacMechAssembler==NULL)
00339     {
00340         Initialise();
00341     }
00342 
00343     BoundaryConditionsContainer<DIM,DIM,1> bcc;
00344     bcc.DefineZeroNeumannOnMeshBoundary(mpElectricsMesh, 0);
00345     mpMonodomainProblem->SetBoundaryConditionsContainer(&bcc);
00346 
00347     // get an electrics assembler from the problem. Note that we don't call
00348     // Solve() on the CardiacProblem class, we do the looping here.
00349     AbstractDynamicAssemblerMixin<DIM,DIM,1>* p_electrics_assembler
00350        = mpMonodomainProblem->CreateAssembler();
00351 
00352     // set up initial voltage etc
00353     Vec voltage=NULL; //This will be set and used later
00354     Vec initial_voltage = mpMonodomainProblem->CreateInitialCondition();
00355 
00356     unsigned num_quad_points = mpCardiacMechAssembler->GetTotalNumQuadPoints();
00357     std::vector<NhsCellularMechanicsOdeSystem> cellmech_systems;
00358     std::vector<double> intracellular_Ca(num_quad_points, 0.0);
00359 
00360     // write the initial position
00361     unsigned counter = 0;
00362 
00363     TimeStepper stepper(0.0, mEndTime, mMechanicsTimeStep);
00364 
00365     unsigned mech_writer_counter = 0;
00366     if (mWriteOutput)
00367     {
00368         mpCardiacMechAssembler->SetWriteOutput();
00369         mpCardiacMechAssembler->WriteOutput(mech_writer_counter);
00370 
00371         if(!mNoElectricsOutput)
00372         {
00373             mpMonodomainProblem->InitialiseWriter();
00374             mpMonodomainProblem->WriteOneStep(stepper.GetTime(), initial_voltage);
00375         }
00376 
00377         if(mIsWatchedLocation)
00378         {
00379             WriteWatchedLocationData(stepper.GetTime(), initial_voltage);
00380         }
00381     }
00382 
00383     while (!stepper.IsTimeAtEnd())
00384     {
00385         LOG(2, "\nCurrent time = " << stepper.GetTime());
00386         std::cout << "\n\n ** Current time = " << stepper.GetTime();
00387 
00388         LOG(2, "  Solving electrics");
00389         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::NON_MECH);
00390         for(unsigned i=0; i<mNumElecTimestepsPerMechTimestep; i++)
00391         {
00392             double current_time = stepper.GetTime() + i*mElectricsTimeStep;
00393             double next_time = stepper.GetTime() + (i+1)*mElectricsTimeStep;
00394 
00395             // solve the electrics
00396             p_electrics_assembler->SetTimes(current_time, next_time, mElectricsTimeStep);
00397             p_electrics_assembler->SetInitialCondition( initial_voltage );
00398 
00399             voltage = p_electrics_assembler->Solve();
00400 
00401             PetscReal min_voltage, max_voltage;
00402             VecMax(voltage,PETSC_NULL,&max_voltage); //the second param is where the index would be returned
00403             VecMin(voltage,PETSC_NULL,&min_voltage);
00404             if(i==0)
00405             {
00406                 LOG(2, "  minimum and maximum voltage is " << min_voltage <<", "<<max_voltage);
00407             }
00408             else if(i==1)
00409             {
00410                 LOG(2, "  ..");
00411             }
00412 
00413             VecDestroy(initial_voltage);
00414             initial_voltage = voltage;
00415         }
00416 
00418 //        p_electrics_assembler->SetMatrixIsNotAssembled();
00419 
00420         // compute Ca_I at each quad point (by interpolation, using the info on which
00421         // electrics element the quad point is in. Then set Ca_I on the mechanics solver
00422         LOG(2, "  Interpolating Ca_I");
00423         for(unsigned i=0; i<mElementAndWeightsForQuadPoints.size(); i++)
00424         {
00425             double interpolated_Ca_I = 0;
00426 
00427             Element<DIM,DIM>& element = *(mpElectricsMesh->GetElement(mElementAndWeightsForQuadPoints[i].ElementNum));
00428             for(unsigned node_index = 0; node_index<element.GetNumNodes(); node_index++)
00429             {
00430                 unsigned global_node_index = element.GetNodeGlobalIndex(node_index);
00431                 double Ca_I_at_node = mpMonodomainProblem->GetPde()->GetCardiacCell(global_node_index)->GetIntracellularCalciumConcentration();
00432                 interpolated_Ca_I += Ca_I_at_node*mElementAndWeightsForQuadPoints[i].Weights(node_index);
00433             }
00434 
00435             intracellular_Ca[i] = interpolated_Ca_I;
00436         }
00437 
00438         LOG(2, "  Setting Ca_I. max value = " << Max(intracellular_Ca));
00439 
00440         // NOTE: HERE WE SHOULD REALLY CHECK WHETHER THE CELL MODELS HAVE Ca_Trop
00441         // AND UPDATE FROM NHS TO CELL_MODEL, BUT NOT SURE HOW TO DO THIS.. (esp for implicit)
00442 
00443         // set [Ca]
00444         mpCardiacMechAssembler->SetIntracellularCalciumConcentrations(intracellular_Ca);
00445         MechanicsEventHandler::EndEvent(MechanicsEventHandler::NON_MECH);
00446 
00447 
00448         // solve the mechanics
00449         LOG(2, "  Solving mechanics ");
00450         //double timestep = std::min(0.01, stepper.GetNextTime()-stepper.GetTime());
00451         mpCardiacMechAssembler->SetWriteOutput(false);
00452 
00453         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ALL_MECH);
00454         mpCardiacMechAssembler->Solve(stepper.GetTime(), stepper.GetNextTime(), mNhsOdeTimeStep);
00455         MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL_MECH);
00456 
00457         LOG(2, "    Number of newton iterations = " << mpCardiacMechAssembler->GetNumNewtonIterations());
00458 
00459         // update the current time
00460         stepper.AdvanceOneTimeStep();
00461         counter++;
00462 
00463         // output the results
00464         MechanicsEventHandler::BeginEvent(MechanicsEventHandler::OUTPUT);
00465         if(mWriteOutput && (counter%WRITE_EVERY_NTH_TIME==0))
00466         {
00467             LOG(2, "  Writing output");
00468             // write deformed position
00469             mech_writer_counter++;
00470             mpCardiacMechAssembler->SetWriteOutput();
00471             mpCardiacMechAssembler->WriteOutput(mech_writer_counter);
00472 
00473             if(!mNoElectricsOutput)
00474             {
00475                 mpMonodomainProblem->mpWriter->AdvanceAlongUnlimitedDimension();
00476                 mpMonodomainProblem->WriteOneStep(stepper.GetTime(), voltage);
00477             }
00478 
00479             if(mIsWatchedLocation)
00480             {
00481                 WriteWatchedLocationData(stepper.GetTime(), voltage);
00482             }
00483         }
00484         MechanicsEventHandler::EndEvent(MechanicsEventHandler::OUTPUT);
00485 
00487 //        // setup the Cinverse data;
00488 //        std::vector<std::vector<double> >& r_c_inverse = NodewiseData<DIM>::Instance()->rGetData();
00489 //        mpCardiacMechAssembler->CalculateCinverseAtNodes(mpElectricsMesh, r_c_inverse);
00490 //
00491 //        // write lambda
00492 //        std::stringstream file_name;
00493 //        file_name << "lambda_" << mech_writer_counter << ".dat";
00494 //        mpCardiacMechAssembler->WriteLambda(mOutputDirectory,file_name.str());
00495 
00496 
00497         // write the total elapsed time..
00498         LogFile::Instance()->WriteElapsedTime("  ");
00499     }
00500 
00501     if ((mWriteOutput) && (!mNoElectricsOutput))
00502     {
00503         HeartConfig::Instance()->Reset();
00504         mpMonodomainProblem->mpWriter->Close();
00505         delete mpMonodomainProblem->mpWriter;
00506 
00507 
00508         // Convert simulation data to Meshalyzer format
00509         //
00510         std::string input_dir = mOutputDirectory+"/electrics";
00511         std::string config_directory = HeartConfig::Instance()->GetOutputDirectory();
00512         HeartConfig::Instance()->SetOutputDirectory(input_dir);
00513         Hdf5ToMeshalyzerConverter converter(input_dir, "voltage");
00514         
00515         // Write mesh in a suitable form for meshalyzer
00516         if (PetscTools::AmMaster())
00517         {
00518             std::string output_directory =  mOutputDirectory + "/electrics/output";
00519             // Write the mesh
00520             MeshalyzerMeshWriter<DIM,DIM> mesh_writer(output_directory, "mesh", false);
00521             mesh_writer.WriteFilesUsingMesh(*mpElectricsMesh);
00522 
00523             // Write the parameters out
00524             HeartConfig::Instance()->Write();
00525         }
00526         
00527         // reset to the default value
00528         HeartConfig::Instance()->SetOutputDirectory(config_directory);
00529         
00530     }
00531 
00532     VecDestroy(voltage);
00533     delete p_electrics_assembler;
00534 
00535     MechanicsEventHandler::EndEvent(MechanicsEventHandler::ALL);
00536 }
00537 
00538 
00539 
00540 template<unsigned DIM>
00541 double CardiacElectroMechanicsProblem<DIM>::Max(std::vector<double>& vec)
00542 {
00543     double max = -1e200;
00544     for(unsigned i=0; i<vec.size(); i++)
00545     {
00546         if(vec[i]>max) max=vec[i];
00547     }
00548     return max;
00549 }
00550 
00551 template<unsigned DIM>
00552 void CardiacElectroMechanicsProblem<DIM>::SetNoElectricsOutput()
00553 {
00554     mNoElectricsOutput = true;
00555 }
00556 
00557 template<unsigned DIM>
00558 void CardiacElectroMechanicsProblem<DIM>::SetWatchedPosition(c_vector<double,DIM> watchedLocation)
00559 {
00560     mIsWatchedLocation = true;
00561     mWatchedLocation = watchedLocation;
00562 }
00563 
00564 template<unsigned DIM>
00565 std::vector<c_vector<double,DIM> >& CardiacElectroMechanicsProblem<DIM>::rGetDeformedPosition()
00566 {
00567     return mpCardiacMechAssembler->rGetDeformedPosition();
00568 }
00569 
00570 
00572 // Explicit instantiation
00574 
00575 //note: 1d incompressible material doesn't make sense
00576 template class CardiacElectroMechanicsProblem<2>;
00577 template class CardiacElectroMechanicsProblem<3>;

Generated on Tue Aug 4 16:10:22 2009 for Chaste by  doxygen 1.5.5