Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM > Class Template Reference

#include <CardiacElectroMechanicsProblem.hpp>

+ Inheritance diagram for CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >:
+ Collaboration diagram for CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >:

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 ()
 
- Public Member Functions inherited from AbstractConductivityModifier< DIM, DIM >
virtual ~AbstractConductivityModifier ()
 
c_matrix< double, SPACE_DIM, SPACE_DIM > & rGetModifiedConductivityTensor (unsigned elementIndex, const c_matrix< double, SPACE_DIM, SPACE_DIM > &rOriginalConductivity, unsigned domainIndex)
 
virtual c_matrix< double, SPACE_DIM, SPACE_DIM > & rCalculateModifiedConductivityTensor (unsigned elementIndex, const c_matrix< double, SPACE_DIM, SPACE_DIM > &rOriginalConductivity, unsigned domainIndex)=0
 

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

  • guess solution
  • this gives the fibre stretch and stretch rate to be set on contraction models
  • integrate contraction models (implicitly if NHS) to get for active tension
  • use this active tension in computing the stress for that guess of the deformation end

Definition at line 94 of file CardiacElectroMechanicsProblem.hpp.

Constructor & Destructor Documentation

◆ CardiacElectroMechanicsProblem()

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 
)

Constructor.

Parameters
compressibilityTypeShould be either INCOMPRESSIBLE or COMPRESSIBLE
electricsProblemTypethe type of electrics problem (MONODOMAIN or BIDOMAIN)
pElectricsMeshMesh on which to solve electrics (Monodomain)
pMechanicsMeshMesh (2nd order) on which to solve mechanics
pCellFactoryfactory to use to create cells
pProblemDefinitionelectro-mechanics problem definition
outputDirectorythe output directory

Definition at line 256 of file CardiacElectroMechanicsProblem.cpp.

References GenericEventHandler< 7, MechanicsEventHandler >::BeginEvent(), CreateElectricsProblem< DIM, PROBLEM_DIM >::Create(), 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().

◆ ~CardiacElectroMechanicsProblem()

template<unsigned DIM, unsigned ELEC_PROB_DIM>
CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::~CardiacElectroMechanicsProblem ( )
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().

Member Function Documentation

◆ DetermineWatchedNodes()

template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::DetermineWatchedNodes ( )
protected

Determine which node is closest to the watched location

Definition at line 60 of file CardiacElectroMechanicsProblem.cpp.

References NEVER_REACHED, OutputFileHandler::OpenOutputFile(), and UNSIGNED_UNSET.

◆ GetMechanicsSolver()

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.

References CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsSolver.

◆ Initialise()

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

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 EXCEPTION, HeartConfig::GetPdeTimeStep(), LogFile::Instance(), HeartConfig::Instance(), NEVER_REACHED, LogFile::Set(), AbstractTetrahedralMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh(), and LogFile::WriteHeader().

◆ Max()

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

Short helper function

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

Definition at line 946 of file CardiacElectroMechanicsProblem.cpp.

◆ OnEndOfTimeStep()

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

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

Parameters
countertime step

Definition at line 291 of file CardiacElectroMechanicsProblem.hpp.

◆ PrepareForSolve()

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

Called in Solve() before the time loop

Definition at line 285 of file CardiacElectroMechanicsProblem.hpp.

◆ rCalculateModifiedConductivityTensor()

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 
)

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

Parameters
elementIndexIndex of current element
rOriginalConductivityReference to the original (for example, undeformed) conductivity tensor
domainIndexUsed 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().

◆ rGetDeformedPosition()

template<unsigned DIM, unsigned ELEC_PROB_DIM>
std::vector< c_vector< double, DIM > > & CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::rGetDeformedPosition ( )
Returns
the current deformed position of the nodes

Definition at line 980 of file CardiacElectroMechanicsProblem.cpp.

◆ SetNoElectricsOutput()

template<unsigned DIM, unsigned ELEC_PROB_DIM>
void CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::SetNoElectricsOutput ( )

Call to not write out voltages

Definition at line 957 of file CardiacElectroMechanicsProblem.cpp.

◆ SetOutputDeformationGradientsAndStress()

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

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
timestephow often to write this. Must be a multiple of the mechanics timestep.

Definition at line 970 of file CardiacElectroMechanicsProblem.cpp.

References EXCEPTION.

◆ SetWatchedPosition()

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

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
watchedLocationlocation (x,y,z) in space. Watched node is the closest to this point.

Definition at line 963 of file CardiacElectroMechanicsProblem.cpp.

◆ Solve()

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

Solve the electromechanics problem

Definition at line 488 of file CardiacElectroMechanicsProblem.cpp.

References TimeStepper::AdvanceOneTimeStep(), PetscTools::Barrier(), GenericEventHandler< 7, MechanicsEventHandler >::BeginEvent(), RelativeTo::ChasteTestOutput, PetscTools::Destroy(), GenericEventHandler< 7, MechanicsEventHandler >::EndEvent(), TimeStepper::GetNextTime(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), HeartConfig::GetOutputDirectory(), HeartConfig::GetPdeTimeStep(), TimeStepper::GetTime(), CommandLineArguments::Instance(), LogFile::Instance(), HeartConfig::Instance(), TimeStepper::IsTimeAtEnd(), CommandLineArguments::OptionExists(), HeartConfig::Reset(), CmguiMeshWriter< ELEMENT_DIM, SPACE_DIM >::SetAdditionalFieldNames(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetInitialCondition(), AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetMatrixIsNotAssembled(), HeartConfig::SetOutputDirectory(), HeartConfig::SetPrintingTimeStep(), CmguiMeshWriter< ELEMENT_DIM, SPACE_DIM >::SetRegionNames(), 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, CmguiDeformedSolutionsWriter< DIM >::WriteCmguiScript(), CmguiDeformedSolutionsWriter< DIM >::WriteDeformationPositions(), LogFile::WriteElapsedTime(), and CmguiDeformedSolutionsWriter< DIM >::WriteInitialMesh().

◆ WriteWatchedLocationData()

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

Write info (x, y, [z], V) for the watched node.

Parameters
timeTime-step now, to write out
voltageVm vector (this is Monodomain)
Todo:
Improve efficiency of this method?

Definition at line 144 of file CardiacElectroMechanicsProblem.cpp.

Friends And Related Symbol Documentation

◆ TestAbstractContractionCellFactory

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
friend class TestAbstractContractionCellFactory
friend

Definition at line 98 of file CardiacElectroMechanicsProblem.hpp.

◆ TestCardiacElectroMechanicsFurtherFunctionality

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
friend class TestCardiacElectroMechanicsFurtherFunctionality
friend

Definition at line 100 of file CardiacElectroMechanicsProblem.hpp.

◆ TestCardiacElectroMechanicsProblem

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
friend class TestCardiacElectroMechanicsProblem
friend

Definition at line 99 of file CardiacElectroMechanicsProblem.hpp.

◆ TestElectroMechanicsExactSolution

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
friend class TestElectroMechanicsExactSolution
friend

Definition at line 101 of file CardiacElectroMechanicsProblem.hpp.

Member Data Documentation

◆ mCompressibilityType

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.

◆ mDeformationGradientsForEachMechanicsElement

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::vector<c_matrix<double,DIM,DIM> > CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationGradientsForEachMechanicsElement
protected

A vector of deformation gradients (each entry a matrix), one for each element in the mechanics mesh

Definition at line 175 of file CardiacElectroMechanicsProblem.hpp.

◆ mDeformationOutputDirectory

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::string CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mDeformationOutputDirectory
protected

◆ mHasBath

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().

◆ mInterpolatedCalciumConcs

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.

◆ mInterpolatedVoltages

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.

◆ mIsWatchedLocation

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
bool CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mIsWatchedLocation
protected

Whether any location has been set to be watched (lots of output for that location

Definition at line 159 of file CardiacElectroMechanicsProblem.hpp.

◆ mModifiedConductivityTensor

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.

◆ mNoElectricsOutput

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

Whether to not write out voltages

Definition at line 153 of file CardiacElectroMechanicsProblem.hpp.

◆ mNumElecTimestepsPerMechTimestep

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.

◆ mNumTimestepsToOutputDeformationGradientsAndStress

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
unsigned CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mNumTimestepsToOutputDeformationGradientsAndStress
protected

How often to print the deformation gradients and stress to file (if ever)

Definition at line 170 of file CardiacElectroMechanicsProblem.hpp.

◆ mOutputDirectory

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
std::string CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mOutputDirectory
protected

Output directory, relative to TEST_OUTPUT

Definition at line 147 of file CardiacElectroMechanicsProblem.hpp.

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

◆ mpCardiacMechSolver

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.

◆ mpElectricsMesh

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
TetrahedralMesh<DIM,DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsMesh
protected

◆ mpElectricsProblem

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
AbstractCardiacProblem<DIM,DIM,ELEC_PROB_DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpElectricsProblem
protected

◆ mpMechanicsMesh

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
QuadraticMesh<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsMesh
protected

◆ mpMechanicsSolver

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
AbstractNonlinearElasticitySolver<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMechanicsSolver
protected

The mechanics solver - a pointer to the part that sees the solid mechanics solver (Object pointed to is the same as with mpCardiacMechSolver)

Definition at line 115 of file CardiacElectroMechanicsProblem.hpp.

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

◆ mpMeshPair

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
FineCoarseMeshPair<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpMeshPair
protected

Class wrapping both meshes, useful for transferring information

Definition at line 144 of file CardiacElectroMechanicsProblem.hpp.

◆ mpProblemDefinition

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
ElectroMechanicsProblemDefinition<DIM>* CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpProblemDefinition
protected

Object containing information about the problem to be solved

Definition at line 138 of file CardiacElectroMechanicsProblem.hpp.

Referenced by CardiacElectroMechProbRegularGeom< DIM >::CardiacElectroMechProbRegularGeom().

◆ mpWatchedLocationFile

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
out_stream CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mpWatchedLocationFile
protected

File where watched location info is written

Definition at line 167 of file CardiacElectroMechanicsProblem.hpp.

◆ mStretchesForEachMechanicsElement

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.

◆ mWatchedElectricsNodeIndex

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
unsigned CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWatchedElectricsNodeIndex
protected

The node in the electrics mesh corresponding to the watched location

Definition at line 163 of file CardiacElectroMechanicsProblem.hpp.

◆ mWatchedLocation

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
c_vector<double,DIM> CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWatchedLocation
protected

The watched location if there is one

Definition at line 161 of file CardiacElectroMechanicsProblem.hpp.

◆ mWatchedMechanicsNodeIndex

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
unsigned CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWatchedMechanicsNodeIndex
protected

The node in the mechanics mesh corresponding to the watched location

Definition at line 165 of file CardiacElectroMechanicsProblem.hpp.

◆ mWriteOutput

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
bool CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::mWriteOutput
protected

◆ WRITE_EVERY_NTH_TIME

template<unsigned DIM, unsigned ELEC_PROB_DIM = 1>
const int CardiacElectroMechanicsProblem< DIM, ELEC_PROB_DIM >::WRITE_EVERY_NTH_TIME = 1
staticprotected

when to write output

Definition at line 156 of file CardiacElectroMechanicsProblem.hpp.


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