Chaste Commit::ca8ccdedf819b6e02855bc0e8e6f50bdecbc5208
ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM > Class Template Reference

#include <ImplicitCardiacMechanicsSolver.hpp>

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

Public Member Functions

 ImplicitCardiacMechanicsSolver (QuadraticMesh< DIM > &rQuadMesh, ElectroMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory)
 
virtual ~ImplicitCardiacMechanicsSolver ()
 
void Solve (double time, double nextTime, double odeTimestep)
 
- Public Member Functions inherited from AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >
 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)
 
void ComputeDeformationGradientAndStretchInEachElement (std::vector< c_matrix< double, DIM, DIM > > &rDeformationGradients, std::vector< double > &rStretches)
 
- Public Member Functions inherited from AbstractCardiacMechanicsSolverInterface< DIM >
 AbstractCardiacMechanicsSolverInterface ()
 
virtual ~AbstractCardiacMechanicsSolverInterface ()
 

Private Member Functions

bool IsImplicitSolver ()
 
void GetActiveTensionAndTensionDerivs (double currentFibreStretch, unsigned currentQuadPointGlobalIndex, bool assembleJacobian, double &rActiveTension, double &rDerivActiveTensionWrtLambda, double &rDerivActiveTensionWrtDLambdaDt)
 

Friends

class TestImplicitCardiacMechanicsSolver
 
class TestExplicitCardiacMechanicsSolver
 

Additional Inherited Members

- Protected Member Functions inherited from AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >
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 ()
 
- Protected Attributes inherited from AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >
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 inherited from AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >
static const unsigned NUM_VERTICES_PER_ELEMENT = ELASTICITY_SOLVER::NUM_VERTICES_PER_ELEMENT
 

Detailed Description

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

Implicit Cardiac Mechanics Solver

The first template parameter should be either IncompressibleNonlinearElasticitySolver or CompressibleNonlinearElasticityAssembler; this will be the class that this class ultimately inherits from.

Solves cardiac mechanics implicitly (together with the contraction models for determining the active tension), taking in the intracellular calcium concentration. See CardiacElectroMechanicsProblem documentation for more detail.

Definition at line 59 of file ImplicitCardiacMechanicsSolver.hpp.

Constructor & Destructor Documentation

◆ ImplicitCardiacMechanicsSolver()

template<class ELASTICITY_SOLVER , unsigned DIM>
ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ImplicitCardiacMechanicsSolver ( 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 39 of file ImplicitCardiacMechanicsSolver.cpp.

◆ ~ImplicitCardiacMechanicsSolver()

Member Function Documentation

◆ GetActiveTensionAndTensionDerivs()

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

A method called by AbstractCardiacMechanicsSolver::AssembleOnElement(), providing the active tension (and other info) at a particular quadrature point. This version uses C to determine the current stretch and stretch rate, and integrates the contraction model ODEs to determine the active tension, and derivatives of the active tension with respect to stretch and stretch rate.

Parameters
currentFibreStretchThe stretch in the fibre direction
currentQuadPointGlobalIndexQuadrature point integrand 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
rDerivActiveTensionWrtDLambdaDtThe returned dT_dLamDot, derivative of active tension wrt stretch rate

Implements AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >.

Definition at line 89 of file ImplicitCardiacMechanicsSolver.cpp.

References DataAtQuadraturePoint_::ContractionModel, EXCEPTION, AbstractContractionModel::GetNextActiveTension(), NEVER_REACHED, AbstractContractionModel::RunDoNotUpdate(), AbstractContractionModel::SetStretchAndStretchRate(), DataAtQuadraturePoint_::Stretch, and DataAtQuadraturePoint_::StretchLastTimeStep.

◆ IsImplicitSolver()

template<class ELASTICITY_SOLVER , unsigned DIM>
bool ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::IsImplicitSolver ( )
inlineprivatevirtual
Returns
true if this solver is an implicit solver (overloaded pure method)

Implements AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >.

Definition at line 65 of file ImplicitCardiacMechanicsSolver.hpp.

◆ Solve()

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

Solve for the deformation using quasi-static nonlinear elasticity. (not dynamic nonlinear elasticity, despite the times taken in - just ONE deformation is solved for. The cell models are integrated implicitly over the time range using the ODE timestep provided, as part of the solve, and updated at the end once the solution has been found, as is lambda.

Parameters
timethe current time
nextTimethe next time
odeTimestepthe ODE timestep

Implements AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >.

Definition at line 57 of file ImplicitCardiacMechanicsSolver.cpp.

References AbstractContractionModel::UpdateStateVariables().

Friends And Related Symbol Documentation

◆ TestExplicitCardiacMechanicsSolver

template<class ELASTICITY_SOLVER , unsigned DIM>
friend class TestExplicitCardiacMechanicsSolver
friend

Definition at line 62 of file ImplicitCardiacMechanicsSolver.hpp.

◆ TestImplicitCardiacMechanicsSolver

template<class ELASTICITY_SOLVER , unsigned DIM>
friend class TestImplicitCardiacMechanicsSolver
friend

Definition at line 61 of file ImplicitCardiacMechanicsSolver.hpp.


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