CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM > Class Template Reference

#include <CardiacElectroMechanicsProblem.hpp>

Inherits AbstractConductivityModifier< DIM, DIM >.

Collaboration diagram for CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 CardiacElectroMechanicsProblem (CompressibilityType compressibilityType, ElectricsProblemType electricsProblemType, TetrahedralMesh< DIM, DIM > *pElectricsMesh, QuadraticMesh< DIM > *pMechanicsMesh, AbstractCardiacCellFactory< DIM > *pCellFactory, ElectroMechanicsProblemDefinition< DIM > *pProblemDefinition, std::string outputDirectory)
virtual ~CardiacElectroMechanicsProblem ()
void Initialise ()
void Solve ()
double Max (std::vector< double > &vec)
void SetNoElectricsOutput ()
void SetWatchedPosition (c_vector< double, DIM > watchedLocation)
void SetOutputDeformationGradientsAndStress (double timestep)
std::vector< c_vector< double,
DIM > > & 
rGetDeformedPosition ()
c_matrix< double, DIM, DIM > & rCalculateModifiedConductivityTensor (unsigned elementIndex, const c_matrix< double, DIM, DIM > &rOriginalConductivity, unsigned domainIndex)
virtual void PrepareForSolve ()
virtual void OnEndOfTimeStep (unsigned counter)
AbstractNonlinearElasticitySolver
< DIM > * 
GetMechanicsSolver ()

Protected Member Functions

void DetermineWatchedNodes ()
void WriteWatchedLocationData (double time, Vec voltage)

Protected Attributes

CompressibilityType mCompressibilityType
AbstractCardiacProblem< DIM,
DIM, ELEC_PROB_DIM > * 
mpElectricsProblem
AbstractCardiacMechanicsSolverInterface
< DIM > * 
mpCardiacMechSolver
AbstractNonlinearElasticitySolver
< DIM > * 
mpMechanicsSolver
unsigned mNumElecTimestepsPerMechTimestep
std::vector< doublemInterpolatedCalciumConcs
std::vector< doublemInterpolatedVoltages
TetrahedralMesh< DIM, DIM > * mpElectricsMesh
QuadraticMesh< DIM > * mpMechanicsMesh
ElectroMechanicsProblemDefinition
< DIM > * 
mpProblemDefinition
bool mHasBath
FineCoarseMeshPair< DIM > * mpMeshPair
std::string mOutputDirectory
std::string mDeformationOutputDirectory
bool mWriteOutput
bool mNoElectricsOutput
bool mIsWatchedLocation
c_vector< double, DIM > mWatchedLocation
unsigned mWatchedElectricsNodeIndex
unsigned mWatchedMechanicsNodeIndex
out_stream mpWatchedLocationFile
unsigned mNumTimestepsToOutputDeformationGradientsAndStress
std::vector< doublemStretchesForEachMechanicsElement
std::vector< c_matrix< double,
DIM, DIM > > 
mDeformationGradientsForEachMechanicsElement
c_matrix< double, DIM, DIM > mModifiedConductivityTensor

Static Protected Attributes

static const int WRITE_EVERY_NTH_TIME = 1

Friends

class TestAbstractContractionCellFactory
class TestCardiacElectroMechanicsProblem
class TestCardiacElectroMechanicsFurtherFunctionality
class TestElectroMechanicsExactSolution

Detailed Description

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
class CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >

CardiacElectroMechanicsProblem

For solving full electro-mechanical problems.

Solves a monodomain problem (diffusion plus cell models) on a (fine) electrics mesh, and a mechanics problem (finite elasticity plus contraction model) on a coarse mesh. An implicit scheme (Jon Whiteley's algorithm) be be used.

For solving problems on regular grids use CardiacElectroMechProbRegularGeom

The implicit algorithm:

Store the position in the electrics mesh of each quad point in the mechanics mesh For every time: Solve the monodomain problem (ie integrate ODEs, solve PDE) Get intracellular [Ca] at each electrics node and interpolate on each mechanics quad point Set [Ca] on each contraction model (one for each point) Solve static finite elasticity problem implicity

Definition at line 94 of file CardiacElectroMechanicsProblem.hpp.


Constructor & Destructor Documentation

template<unsigned DIM, unsigned ELEC_PROB_DIM>
CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::CardiacElectroMechanicsProblem ( CompressibilityType  compressibilityType,
ElectricsProblemType  electricsProblemType,
TetrahedralMesh< DIM, DIM > *  pElectricsMesh,
QuadraticMesh< DIM > *  pMechanicsMesh,
AbstractCardiacCellFactory< DIM > *  pCellFactory,
ElectroMechanicsProblemDefinition< DIM > *  pProblemDefinition,
std::string  outputDirectory 
) [inline]

Constructor.

Parameters:
compressibilityType Should be either INCOMPRESSIBLE or COMPRESSIBLE
electricsProblemType the type of electrics problem (MONODOMAIN or BIDOMAIN)
pElectricsMesh Mesh on which to solve electrics (Monodomain)
pMechanicsMesh Mesh (2nd order) on which to solve mechanics
pCellFactory factory to use to create cells
pProblemDefinition electro-mechanics problem definition
outputDirectory the output directory

Definition at line 256 of file CardiacElectroMechanicsProblem.cpp.

References GenericEventHandler< 7, MechanicsEventHandler >::BeginEvent(), GenericEventHandler< 16, HeartEventHandler >::Disable(), HeartConfig::Instance(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationOutputDirectory, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mHasBath, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mOutputDirectory, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsProblem, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWriteOutput, GenericEventHandler< 7, MechanicsEventHandler >::Reset(), HeartConfig::SetOutputDirectory(), and HeartConfig::SetOutputFilenamePrefix().

template<unsigned DIM, unsigned ELEC_PROB_DIM>
CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::~CardiacElectroMechanicsProblem (  )  [inline, virtual]

Delete allocated memory and close the watched location file

NOTE if SetWatchedLocation but not Initialise has been called, mpWatchedLocationFile will be uninitialised and using it will cause a seg fault. Hence the mpMechanicsMesh!=NULL it is true if Initialise has been called.

Definition at line 327 of file CardiacElectroMechanicsProblem.cpp.

References LogFile::Close(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpCardiacMechSolver, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsProblem, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsMesh, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMeshPair, and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpWatchedLocationFile.


Member Function Documentation

template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::DetermineWatchedNodes (  )  [inline, protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
AbstractNonlinearElasticitySolver<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::GetMechanicsSolver (  )  [inline]

Get the mechanics solver used in the solve. Only call after a solve.

Returns:
The mechanics (nonlinear elasticity) PDE solver.

Definition at line 300 of file CardiacElectroMechanicsProblem.hpp.

template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Initialise (  )  [inline]

Initialise the class. Initialises the MonodomainProblem and sets up the electrics mesh to mechanics mesh data.

Todo:
This is fragile: check how the TimeStepper does it, and possibly refactor the behaviour there into a static helper method if it isn't already.

Definition at line 345 of file CardiacElectroMechanicsProblem.cpp.

References CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::DetermineWatchedNodes(), EXCEPTION, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), HeartConfig::GetPdeTimeStep(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::GetTissue(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Initialise(), LogFile::Instance(), HeartConfig::Instance(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mCompressibilityType, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationGradientsForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationOutputDirectory, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mInterpolatedCalciumConcs, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mInterpolatedVoltages, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNumElecTimestepsPerMechTimestep, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mOutputDirectory, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpCardiacMechSolver, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsMesh, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsProblem, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsMesh, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsSolver, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMeshPair, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpProblemDefinition, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mStretchesForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWriteOutput, NEVER_REACHED, LogFile::Set(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetMesh(), AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh(), and LogFile::WriteHeader().

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM>
double CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Max ( std::vector< double > &  vec  )  [inline]

Short helper function

Returns:
the max of a std::vector
Parameters:
vec a vector of doubles

Definition at line 950 of file CardiacElectroMechanicsProblem.cpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
virtual void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::OnEndOfTimeStep ( unsigned  counter  )  [inline, virtual]

Called in Solve() at the end of every time step

Parameters:
counter time step

Definition at line 291 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
virtual void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::PrepareForSolve (  )  [inline, virtual]

Called in Solve() before the time loop

Definition at line 285 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM>
c_matrix< double, DIM, DIM > & CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::rCalculateModifiedConductivityTensor ( unsigned  elementIndex,
const c_matrix< double, DIM, DIM > &  rOriginalConductivity,
unsigned  domainIndex 
) [inline]

The implementation of the pure method defined in the base class AbstractConductivityModifier. The tissue class will call this method.

Parameters:
elementIndex Index of current element
rOriginalConductivity Reference to the original (for example, undeformed) conductivity tensor
domainIndex Used to tailor modification to the domain. 0 = intracellular, 1 = extracellular, 2 = second intracellular (tridomain)
Returns:
Reference to a modified conductivity tensor.

Definition at line 169 of file CardiacElectroMechanicsProblem.cpp.

References Inverse(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationGradientsForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mModifiedConductivityTensor, and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMeshPair.

template<unsigned DIM, unsigned ELEC_PROB_DIM>
std::vector< c_vector< double, DIM > > & CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::rGetDeformedPosition (  )  [inline]
template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::SetNoElectricsOutput (  )  [inline]

Call to not write out voltages

Definition at line 961 of file CardiacElectroMechanicsProblem.cpp.

References CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNoElectricsOutput.

template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::SetOutputDeformationGradientsAndStress ( double  timestep  )  [inline]

Call this for a files containing the deformation gradients (F), evaluated at the centroids of element, and the 2nd PK stress for each element (averaged over the values at the quadrature points of that element), to be written,

Parameters:
timestep how often to write this. Must be a multiple of the mechanics timestep.

Definition at line 974 of file CardiacElectroMechanicsProblem.cpp.

References EXCEPTION, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNumTimestepsToOutputDeformationGradientsAndStress, and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpProblemDefinition.

template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::SetWatchedPosition ( c_vector< double, DIM >  watchedLocation  )  [inline]

Set a location to be watched - for which lots of output is given. Should correspond to nodes in both meshes.

The watched file will have rows that look like: time x_pos y_pos [z_pos] voltage Ca_i_conc.

NOTE: for the Calcium - assumes LUO_RUDY IS USED

Parameters:
watchedLocation location (x,y,z) in space. Watched node is the closest to this point.

Definition at line 967 of file CardiacElectroMechanicsProblem.cpp.

References CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mIsWatchedLocation, and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWatchedLocation.

template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve (  )  [inline]

Solve the electromechanics problem

Definition at line 488 of file CardiacElectroMechanicsProblem.cpp.

References Hdf5DataWriter::AdvanceAlongUnlimitedDimension(), PetscTools::Barrier(), GenericEventHandler< 7, MechanicsEventHandler >::BeginEvent(), RelativeTo::ChasteTestOutput, Hdf5DataWriter::Close(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::CreateInitialCondition(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::CreateSolver(), DistributedVectorFactory::CreateVec(), PetscTools::Destroy(), GenericEventHandler< 7, MechanicsEventHandler >::EndEvent(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetDistributedVectorFactory(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetElement(), DistributedVectorFactory::GetHigh(), DistributedVectorFactory::GetLow(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), HeartConfig::GetOutputDirectory(), HeartConfig::GetPdeTimeStep(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::GetTissue(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Initialise(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::InitialiseWriter(), LogFile::Instance(), HeartConfig::Instance(), CommandLineArguments::Instance(), DistributedVectorFactory::IsGlobalIndexLocal(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Max(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationGradientsForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mHasBath, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mInterpolatedCalciumConcs, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mInterpolatedVoltages, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mIsWatchedLocation, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNoElectricsOutput, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNumElecTimestepsPerMechTimestep, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNumTimestepsToOutputDeformationGradientsAndStress, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mOutputDirectory, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpCardiacMechSolver, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsMesh, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsProblem, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsMesh, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsSolver, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMeshPair, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpProblemDefinition, AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpWriter, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mStretchesForEachMechanicsElement, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWriteOutput, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::OnEndOfTimeStep(), CommandLineArguments::OptionExists(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::PrepareForSolve(), HeartConfig::Reset(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::rGetDeformedPosition(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetBoundaryConditionsContainer(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetInitialCondition(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetMatrixIsNotAssembled(), HeartConfig::SetOutputDirectory(), HeartConfig::SetPrintingTimeStep(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetTimes(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetTimeStep(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve(), UNSIGNED_UNSET, CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::WRITE_EVERY_NTH_TIME, LogFile::WriteElapsedTime(), AbstractCardiacProblem< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::WriteOneStep(), and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::WriteWatchedLocationData().

template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::WriteWatchedLocationData ( double  time,
Vec  voltage 
) [inline, protected]

Member Data Documentation

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
CompressibilityType CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mCompressibilityType [protected]

Either COMPRESSIBLE or INCOMPRESSIBLE

Definition at line 105 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Initialise().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::vector<c_matrix<double,DIM,DIM> > CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationGradientsForEachMechanicsElement [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::string CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationOutputDirectory [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
bool CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mHasBath [protected]

Whether the mesh has a bath (non-active) region or not. False by default.

Definition at line 141 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::CardiacElectroMechanicsProblem(), and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::vector<double> CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mInterpolatedCalciumConcs [protected]

A cache for the interpolated calcium concentrations from electrics to mechanics mesh. Memory is allocated within Initialise(). Filled in during Solve() and passed on to the mechanics solver

Definition at line 124 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::vector<double> CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mInterpolatedVoltages [protected]

A cache for the interpolated voltages from electrics to mechanics mesh. Memory is allocated within Initialise(). Filled in during Solve() and passed on to the mechanics solver

Definition at line 130 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
bool CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mIsWatchedLocation [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
c_matrix<double,DIM,DIM> CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mModifiedConductivityTensor [protected]

Somewhere to store the modified conductivity tensor

Definition at line 178 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::rCalculateModifiedConductivityTensor().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
bool CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNoElectricsOutput [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
unsigned CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNumElecTimestepsPerMechTimestep [protected]

The number of electrics timesteps per mechanics timestep

Definition at line 118 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
unsigned CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNumTimestepsToOutputDeformationGradientsAndStress [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::string CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mOutputDirectory [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
AbstractCardiacMechanicsSolverInterface<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpCardiacMechSolver [protected]

The mechanics solver - a pointer to the part that sees the cardiac mechanics interface bit. (Object pointed to is the same as with mpMechanicsSolver)

Definition at line 112 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Initialise(), CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve(), and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::~CardiacElectroMechanicsProblem().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
TetrahedralMesh<DIM,DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsMesh [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
AbstractCardiacProblem<DIM,DIM,ELEC_PROB_DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsProblem [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
QuadraticMesh<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsMesh [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
AbstractNonlinearElasticitySolver<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsSolver [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
FineCoarseMeshPair<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMeshPair [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
ElectroMechanicsProblemDefinition<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpProblemDefinition [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
out_stream CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpWatchedLocationFile [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::vector<double> CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mStretchesForEachMechanicsElement [protected]

A vector of stretches (in the fibre direction), one for each element in the mechanics mesh

Definition at line 173 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Initialise(), and CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
unsigned CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWatchedElectricsNodeIndex [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
c_vector<double,DIM> CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWatchedLocation [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
unsigned CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWatchedMechanicsNodeIndex [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
bool CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWriteOutput [protected]
template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
const int CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::WRITE_EVERY_NTH_TIME = 1 [static, protected]

when to write output

Definition at line 156 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::Solve().


The documentation for this class was generated from the following files:

Generated by  doxygen 1.6.2