Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM > Class Template Referenceabstract

#include <AbstractCardiacMechanicsSolver.hpp>

+ Inheritance diagram for AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >:
+ Collaboration diagram for AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >:

Public Member Functions

 AbstractCardiacMechanicsSolver (QuadraticMesh< DIM > &rQuadMesh, ElectroMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory)
 
virtual ~AbstractCardiacMechanicsSolver ()
 
void SetFineCoarseMeshPair (FineCoarseMeshPair< DIM > *pMeshPair)
 
unsigned GetTotalNumQuadPoints ()
 
virtual GaussianQuadratureRule< DIM > * GetQuadratureRule ()
 
std::map< unsigned, DataAtQuadraturePoint > & rGetQuadPointToDataAtQuadPointMap ()
 
void SetConstantFibreSheetDirections (const c_matrix< double, DIM, DIM > &rFibreSheetMatrix)
 
void SetVariableFibreSheetDirections (const FileFinder &rOrthoFile, bool definedPerQuadraturePoint)
 
void SetCalciumAndVoltage (std::vector< double > &rCalciumConcentrations, std::vector< double > &rVoltages)
 
virtual void Solve (double time, double nextTime, double odeTimestep)=0
 
void ComputeDeformationGradientAndStretchInEachElement (std::vector< c_matrix< double, DIM, DIM > > &rDeformationGradients, std::vector< double > &rStretches)
 
- Public Member Functions inherited from AbstractCardiacMechanicsSolverInterface< DIM >
 AbstractCardiacMechanicsSolverInterface ()
 
virtual ~AbstractCardiacMechanicsSolverInterface ()
 

Protected Member Functions

virtual bool IsImplicitSolver ()=0
 
void AddActiveStressAndStressDerivative (c_matrix< double, DIM, DIM > &rC, unsigned elementIndex, unsigned currentQuadPointGlobalIndex, c_matrix< double, DIM, DIM > &rT, FourthOrderTensor< DIM, DIM, DIM, DIM > &rDTdE, bool addToDTdE)
 
void SetupChangeOfBasisMatrix (unsigned elementIndex, unsigned currentQuadPointGlobalIndex)
 
void Initialise ()
 
virtual void GetActiveTensionAndTensionDerivs (double currentFibreStretch, unsigned currentQuadPointGlobalIndex, bool assembleJacobian, double &rActiveTension, double &rDerivActiveTensionWrtLambda, double &rDerivActiveTensionWrtDLambdaDt)=0
 

Protected Attributes

std::map< unsigned, DataAtQuadraturePointmQuadPointToDataAtQuadPointMap
 
std::map< unsigned, DataAtQuadraturePoint >::iterator mMapIterator
 
FineCoarseMeshPair< DIM > * mpMeshPair
 
unsigned mTotalQuadPoints
 
double mCurrentTime
 
double mNextTime
 
double mOdeTimestep
 
c_matrix< double, DIM, DIM > mConstantFibreSheetDirections
 
std::vector< c_matrix< double, DIM, DIM > > * mpVariableFibreSheetDirections
 
bool mFibreSheetDirectionsDefinedByQuadraturePoint
 
c_vector< double, DIM > mCurrentElementFibreDirection
 
c_vector< double, DIM > mCurrentElementSheetDirection
 
c_vector< double, DIM > mCurrentElementSheetNormalDirection
 
ElectroMechanicsProblemDefinition< DIM > & mrElectroMechanicsProblemDefinition
 

Static Protected Attributes

static const unsigned NUM_VERTICES_PER_ELEMENT = ELASTICITY_SOLVER::NUM_VERTICES_PER_ELEMENT
 

Detailed Description

template<class ELASTICITY_SOLVER, unsigned DIM>
class AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >

AbstractCardiacMechanicsSolver

Base class to implicit and explicit cardiac mechanics solvers. Inherits from IncompressibleNonlinearElasticitySolver or CompressibleNonlinearElasticityAssembler (depending on what the template parameter ELASTICITY_SOLVER is), and also from AbstractCardiacMechanicsSolverInterface which just declares this classes main public methods.

Overloads AddActiveStressAndStressDerivative() which adds on the active tension term to the stress. The child classes hold the contraction models and need to implement a method for getting the active tension from the model.

Definition at line 81 of file AbstractCardiacMechanicsSolver.hpp.

Constructor & Destructor Documentation

◆ AbstractCardiacMechanicsSolver()

template<class ELASTICITY_SOLVER , unsigned DIM>
AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::AbstractCardiacMechanicsSolver ( QuadraticMesh< DIM > &  rQuadMesh,
ElectroMechanicsProblemDefinition< DIM > &  rProblemDefinition,
std::string  outputDirectory 
)

Constructor

Parameters
rQuadMeshA reference to the mesh.
rProblemDefinitionObject defining body force and boundary conditions
outputDirectoryThe output directory, relative to TEST_OUTPUT

Definition at line 41 of file AbstractCardiacMechanicsSolver.cpp.

◆ ~AbstractCardiacMechanicsSolver()

Member Function Documentation

◆ AddActiveStressAndStressDerivative()

template<class ELASTICITY_SOLVER , unsigned DIM>
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::AddActiveStressAndStressDerivative ( c_matrix< double, DIM, DIM > &  rC,
unsigned  elementIndex,
unsigned  currentQuadPointGlobalIndex,
c_matrix< double, DIM, DIM > &  rT,
FourthOrderTensor< DIM, DIM, DIM, DIM > &  rDTdE,
bool  addToDTdE 
)
protected

Overloaded AddActiveStressAndStressDerivative(), which calls on the contraction model to get the active stress and add it on to the stress tensor

Parameters
rCThe Lagrangian deformation tensor (F^T F)
elementIndexIndex of the current element
currentQuadPointGlobalIndexThe index (assuming an outer loop over elements and an inner loop over quadrature points), of the current quadrature point.
rTThe stress to be added to
rDTdEthe stress derivative to be added to, assuming the final parameter is true
addToDTdEA boolean flag saying whether the stress derivative is required or not.

Definition at line 208 of file AbstractCardiacMechanicsSolver.cpp.

References Determinant(), and Inverse().

◆ ComputeDeformationGradientAndStretchInEachElement()

template<class ELASTICITY_SOLVER , unsigned DIM>
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement ( std::vector< c_matrix< double, DIM, DIM > > &  rDeformationGradients,
std::vector< double > &  rStretches 
)
virtual

Compute the deformation gradient, and stretch in the fibre direction, for each element in the mesh. Note: using quadratic interpolation for position, the deformation gradients and stretches actually vary linearly in each element. However, for computational efficiency reasons, when computing deformation gradients and stretches to pass back to the electrophysiology solver, we just assume they are constant in each element (ie ignoring the quadratic correction to the displacement). This means that the (const) deformation gradient and stretch for each element can be computed in advance and stored, and we don't have to worry about interpolation onto the precise location of the cell-model (electrics-mesh) node, just which element it is in, and ditto the electric mesh element centroid.

To compute this (elementwise-)constant F (and from it the constant stretch), we just have to compute F using the deformed positions at the vertices only, with linear bases, rather than all the nodes and quadratic bases.

Parameters
rDeformationGradientsA reference of a std::vector in which the deformation gradient in each element will be returned. Must be allocated prior to being passed in.
rStretchesA reference of a std::vector in which the stretch in each element will be returned. Must be allocated prior to being passed in.

Implements AbstractCardiacMechanicsSolverInterface< DIM >.

Definition at line 418 of file AbstractCardiacMechanicsSolver.cpp.

References LinearBasisFunction< ELEMENT_DIM >::ComputeTransformedBasisFunctionDerivatives(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex().

◆ GetActiveTensionAndTensionDerivs()

template<class ELASTICITY_SOLVER , unsigned DIM>
virtual void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs ( double  currentFibreStretch,
unsigned  currentQuadPointGlobalIndex,
bool  assembleJacobian,
double rActiveTension,
double rDerivActiveTensionWrtLambda,
double rDerivActiveTensionWrtDLambdaDt 
)
protectedpure virtual

Pure method called in AbstractCardiacMechanicsSolver::AddActiveStressAndStressDerivative(), which needs to provide the active tension (and other info if implicit (if the contraction model depends on stretch or stretch rate)) at a particular quadrature point. Takes in the current fibre stretch.

Parameters
currentFibreStretchThe stretch in the fibre direction
currentQuadPointGlobalIndexquadrature point the integrand is currently being evaluated at in AssembleOnElement
assembleJacobianA bool stating whether to assemble the Jacobian matrix.
rActiveTensionThe returned active tension
rDerivActiveTensionWrtLambdaThe returned dT_dLam, derivative of active tension wrt stretch. Only should be set in implicit solvers
rDerivActiveTensionWrtDLambdaDtThe returned dT_dLamDot, derivative of active tension wrt stretch rate. Only should be set in implicit solver

Implemented in ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >, and ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >.

◆ GetQuadratureRule()

template<class ELASTICITY_SOLVER , unsigned DIM>
virtual GaussianQuadratureRule< DIM > * AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetQuadratureRule ( )
inlinevirtual
Returns
the quadrature rule used in the elements.

Implements AbstractCardiacMechanicsSolverInterface< DIM >.

Definition at line 247 of file AbstractCardiacMechanicsSolver.hpp.

◆ GetTotalNumQuadPoints()

template<class ELASTICITY_SOLVER , unsigned DIM>
unsigned AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetTotalNumQuadPoints ( )
inlinevirtual
Returns
the total number of quad points in the mesh. Pure, implemented in concrete solver

Implements AbstractCardiacMechanicsSolverInterface< DIM >.

Definition at line 241 of file AbstractCardiacMechanicsSolver.hpp.

References AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mTotalQuadPoints.

◆ Initialise()

template<class ELASTICITY_SOLVER , unsigned DIM>
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Initialise ( )
protectedvirtual

◆ IsImplicitSolver()

template<class ELASTICITY_SOLVER , unsigned DIM>
virtual bool AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::IsImplicitSolver ( )
protectedpure virtual
Returns
whether the solver is implicit or not (i.e. whether the contraction model depends on lambda (and depends on lambda at the current time)). For whether dTa_dLam dependent terms need to be added to the Jacbobian

Implemented in ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >, and ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >.

◆ rGetQuadPointToDataAtQuadPointMap()

template<class ELASTICITY_SOLVER , unsigned DIM>
std::map< unsigned, DataAtQuadraturePoint > & AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::rGetQuadPointToDataAtQuadPointMap ( )
inline
Returns
access mQuadPointToDataAtQuadPointMap. See doxygen for this variable

Definition at line 255 of file AbstractCardiacMechanicsSolver.hpp.

References AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mQuadPointToDataAtQuadPointMap.

◆ SetCalciumAndVoltage()

template<class ELASTICITY_SOLVER , unsigned DIM>
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetCalciumAndVoltage ( std::vector< double > &  rCalciumConcentrations,
std::vector< double > &  rVoltages 
)
virtual

Set the intracellular Calcium concentrations and voltages at each quad point. Pure.

Implicit solvers (for contraction models which are functions of stretch (and maybe stretch rate) would integrate the contraction model with this Ca/V/t using the current stretch (ie inside AssembleOnElement, ie inside GetActiveTensionAndTensionDerivs). Explicit solvers (for contraction models which are NOT functions of stretch can immediately integrate the contraction models to get the active tension.

Parameters
rCalciumConcentrationsReference to a vector of intracellular calcium concentrations at each quadrature point
rVoltagesReference to a vector of voltages at each quadrature point
Todo:
#1828 / #1211 don't pass in entire vector

Implements AbstractCardiacMechanicsSolverInterface< DIM >.

Definition at line 165 of file AbstractCardiacMechanicsSolver.cpp.

References ContractionModelInputParameters_::intracellularCalciumConcentration, and ContractionModelInputParameters_::voltage.

◆ SetConstantFibreSheetDirections()

template<class ELASTICITY_SOLVER , unsigned DIM>
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetConstantFibreSheetDirections ( const c_matrix< double, DIM, DIM > &  rFibreSheetMatrix)
virtual

Set a constant fibre-sheet-normal direction (a matrix) to something other than the default (fibres in X-direction, sheet in the XY plane)

Parameters
rFibreSheetMatrixThe fibre-sheet-normal matrix (fibre dir the first column, normal-to-fibre-in sheet in second column, sheet-normal in third column).

Implements AbstractCardiacMechanicsSolverInterface< DIM >.

Definition at line 517 of file AbstractCardiacMechanicsSolver.cpp.

References EXCEPTION.

◆ SetFineCoarseMeshPair()

template<class ELASTICITY_SOLVER , unsigned DIM>
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetFineCoarseMeshPair ( FineCoarseMeshPair< DIM > *  pMeshPair)
virtual

Sets the fine-coarse mesh pair object so that the solver knows about electrics too. It checks that the coarse mesh of the fine-mesh pair has the same number of elements as the quad mesh of this object and throws an exception otherwise.

Parameters
pMeshPairthe FineCoarseMeshPair object to be set

Implements AbstractCardiacMechanicsSolverInterface< DIM >.

Definition at line 134 of file AbstractCardiacMechanicsSolver.cpp.

References EXCEPTION, FineCoarseMeshPair< DIM >::GetCoarseMesh(), and AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements().

◆ SetupChangeOfBasisMatrix()

template<class ELASTICITY_SOLVER , unsigned DIM>
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetupChangeOfBasisMatrix ( unsigned  elementIndex,
unsigned  currentQuadPointGlobalIndex 
)
protected

Over-ridden method which sets up an internal variable in the parent class, using the provided fibre-sheet direction information.

Parameters
elementIndexelement global index
currentQuadPointGlobalIndexquad point global index

Definition at line 190 of file AbstractCardiacMechanicsSolver.cpp.

◆ SetVariableFibreSheetDirections()

template<class ELASTICITY_SOLVER , unsigned DIM>
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetVariableFibreSheetDirections ( const FileFinder rOrthoFile,
bool  definedPerQuadraturePoint 
)
virtual

Set a variable fibre-sheet-normal direction (matrices), from file. If the second parameter is false, there should be one fibre-sheet definition for each element; otherwise there should be one fibre-sheet definition for each *quadrature point* in the mesh. In the first case, the file should be a .ortho file (ie each line has the fibre dir, sheet dir, normal dir for that element), in the second it should have .orthoquad as the format.

Parameters
rOrthoFilethe file containing the fibre/sheet directions
definedPerQuadraturePointwhether the fibre-sheet definitions are for each quadrature point in the mesh (if not, one for each element is assumed).

Implements AbstractCardiacMechanicsSolverInterface< DIM >.

Definition at line 487 of file AbstractCardiacMechanicsSolver.cpp.

References EXCEPTION, FileFinder::GetAbsolutePath(), FibreReader< DIM >::GetFibreSheetAndNormalMatrix(), and FibreReader< DIM >::GetNumLinesOfData().

◆ Solve()

template<class ELASTICITY_SOLVER , unsigned DIM>
virtual void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve ( double  time,
double  nextTime,
double  odeTimestep 
)
pure virtual

Solve for the deformation, integrating the contraction model ODEs.

Parameters
timethe current time
nextTimethe next time
odeTimestepthe ODE timestep

Implements AbstractCardiacMechanicsSolverInterface< DIM >.

Implemented in ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >, and ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >.

Member Data Documentation

◆ mConstantFibreSheetDirections

template<class ELASTICITY_SOLVER , unsigned DIM>
c_matrix<double,DIM,DIM> AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mConstantFibreSheetDirections
protected

The fibre-sheet-normal directions (in a matrix), if constant (defaults to the identity, ie fibres in the X-direction, sheet in the XY plane)

Definition at line 122 of file AbstractCardiacMechanicsSolver.hpp.

◆ mCurrentElementFibreDirection

template<class ELASTICITY_SOLVER , unsigned DIM>
c_vector<double,DIM> AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentElementFibreDirection
protected

The fibre direction for the current element being assembled on

Definition at line 137 of file AbstractCardiacMechanicsSolver.hpp.

◆ mCurrentElementSheetDirection

template<class ELASTICITY_SOLVER , unsigned DIM>
c_vector<double,DIM> AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentElementSheetDirection
protected

The sheet direction for the current element being assembled on

Definition at line 140 of file AbstractCardiacMechanicsSolver.hpp.

◆ mCurrentElementSheetNormalDirection

template<class ELASTICITY_SOLVER , unsigned DIM>
c_vector<double,DIM> AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentElementSheetNormalDirection
protected

The sheet normal direction for the current element being assembled on

Definition at line 143 of file AbstractCardiacMechanicsSolver.hpp.

◆ mCurrentTime

template<class ELASTICITY_SOLVER , unsigned DIM>
double AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentTime
protected

Current time

Definition at line 110 of file AbstractCardiacMechanicsSolver.hpp.

◆ mFibreSheetDirectionsDefinedByQuadraturePoint

template<class ELASTICITY_SOLVER , unsigned DIM>
bool AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mFibreSheetDirectionsDefinedByQuadraturePoint
protected

Whether the fibre-sheet directions that where read in where define per element or per quadrature point. Only valid if mpVariableFibreSheetDirections!=NULL

Definition at line 134 of file AbstractCardiacMechanicsSolver.hpp.

◆ mMapIterator

template<class ELASTICITY_SOLVER , unsigned DIM>
std::map<unsigned,DataAtQuadraturePoint>::iterator AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mMapIterator
protected

An iterator to the map, used to avoid having to repeatedly search the map.

Definition at line 101 of file AbstractCardiacMechanicsSolver.hpp.

◆ mNextTime

template<class ELASTICITY_SOLVER , unsigned DIM>
double AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mNextTime
protected

Time to which the solver has been asked to solve to

Definition at line 113 of file AbstractCardiacMechanicsSolver.hpp.

◆ mOdeTimestep

template<class ELASTICITY_SOLVER , unsigned DIM>
double AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mOdeTimestep
protected

Time used to integrate the contraction model

Definition at line 116 of file AbstractCardiacMechanicsSolver.hpp.

◆ mpMeshPair

template<class ELASTICITY_SOLVER , unsigned DIM>
FineCoarseMeshPair<DIM>* AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpMeshPair
protected

A mesh pair object that can be set by the user to inform the solver about the electrics mesh.

Definition at line 104 of file AbstractCardiacMechanicsSolver.hpp.

◆ mpVariableFibreSheetDirections

template<class ELASTICITY_SOLVER , unsigned DIM>
std::vector<c_matrix<double,DIM,DIM> >* AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpVariableFibreSheetDirections
protected

The fibre-sheet-normal directions (matrices), one for each element. Only non-NULL if SetVariableFibreSheetDirections() is called, if not mConstantFibreSheetDirections is used instead

Definition at line 128 of file AbstractCardiacMechanicsSolver.hpp.

◆ mQuadPointToDataAtQuadPointMap

template<class ELASTICITY_SOLVER , unsigned DIM>
std::map<unsigned,DataAtQuadraturePoint> AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mQuadPointToDataAtQuadPointMap
protected

A map from the index of a quadrature point to the data (contraction model, stretch, stretch at the last time-step) at that quad point. Note that there is no vector of all the quadrature points of the mesh; the quad point index is the index that would be obtained by looping over elements and then looping over quad points.

DISTRIBUTED - only holds data for the quad points within elements owned by this process.

Definition at line 96 of file AbstractCardiacMechanicsSolver.hpp.

Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::rGetQuadPointToDataAtQuadPointMap().

◆ mrElectroMechanicsProblemDefinition

template<class ELASTICITY_SOLVER , unsigned DIM>
ElectroMechanicsProblemDefinition<DIM>& AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mrElectroMechanicsProblemDefinition
protected

This class contains all the information about the electro mechanics problem (except the material law)

Definition at line 149 of file AbstractCardiacMechanicsSolver.hpp.

◆ mTotalQuadPoints

template<class ELASTICITY_SOLVER , unsigned DIM>
unsigned AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mTotalQuadPoints
protected

Total number of quad points in the (mechanics) mesh

Definition at line 107 of file AbstractCardiacMechanicsSolver.hpp.

Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetTotalNumQuadPoints().

◆ NUM_VERTICES_PER_ELEMENT

template<class ELASTICITY_SOLVER , unsigned DIM>
const unsigned AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::NUM_VERTICES_PER_ELEMENT = ELASTICITY_SOLVER::NUM_VERTICES_PER_ELEMENT
staticprotected

Useful const from base class

Definition at line 84 of file AbstractCardiacMechanicsSolver.hpp.


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