#include <AbstractCardiacMechanicsSolver.hpp>
Inherits ELASTICITY_SOLVER, and AbstractCardiacMechanicsSolverInterface< DIM >.
Inherited by ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >, and ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >.
Public Member Functions | |
AbstractCardiacMechanicsSolver (QuadraticMesh< DIM > &rQuadMesh, SolidMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory) | |
~AbstractCardiacMechanicsSolver () | |
unsigned | GetTotalNumQuadPoints () |
virtual GaussianQuadratureRule < DIM > * | GetQuadratureRule () |
void | SetConstantFibreSheetDirections (const c_matrix< double, DIM, DIM > &rFibreSheetMatrix) |
void | SetVariableFibreSheetDirections (std::string orthoFile, 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) |
Protected Member Functions | |
virtual bool | IsImplicitSolver ()=0 |
void | ComputeStressAndStressDerivative (AbstractMaterialLaw< DIM > *pMaterialLaw, c_matrix< double, DIM, DIM > &rC, c_matrix< double, DIM, DIM > &rInvC, double pressure, unsigned elementIndex, unsigned currentQuadPointGlobalIndex, c_matrix< double, DIM, DIM > &rT, FourthOrderTensor< DIM, DIM, DIM, DIM > &rDTdE, bool computeDTdE) |
virtual void | GetActiveTensionAndTensionDerivs (double currentFibreStretch, unsigned currentQuadPointGlobalIndex, bool assembleJacobian, double &rActiveTension, double &rDerivActiveTensionWrtLambda, double &rDerivActiveTensionWrtDLambdaDt)=0 |
Protected Attributes | |
std::vector < AbstractContractionModel * > | mContractionModelSystems |
std::vector< double > | mStretches |
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_matrix< double, DIM, DIM > * | mpCurrentElementFibreSheetMatrix |
c_vector< double, DIM > | mCurrentElementFibreDirection |
Static Protected Attributes | |
static const unsigned | NUM_VERTICES_PER_ELEMENT = ELASTICITY_SOLVER::NUM_VERTICES_PER_ELEMENT |
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 ComputeStressAndStressDerivative() 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 54 of file AbstractCardiacMechanicsSolver.hpp.
AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::AbstractCardiacMechanicsSolver | ( | QuadraticMesh< DIM > & | rQuadMesh, | |
SolidMechanicsProblemDefinition< DIM > & | rProblemDefinition, | |||
std::string | outputDirectory | |||
) | [inline] |
Constructor
rQuadMesh | A reference to the mesh. | |
rProblemDefinition | Object defining body force and boundary conditions | |
outputDirectory | The output directory, relative to TEST_OUTPUT |
Definition at line 263 of file AbstractCardiacMechanicsSolver.hpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mConstantFibreSheetDirections, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpVariableFibreSheetDirections, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mStretches, and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mTotalQuadPoints.
AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::~AbstractCardiacMechanicsSolver | ( | ) | [inline] |
Destructor just deletes memory if it was allocated
Definition at line 291 of file AbstractCardiacMechanicsSolver.hpp.
References AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpVariableFibreSheetDirections.
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement | ( | std::vector< c_matrix< double, DIM, DIM > > & | rDeformationGradients, | |
std::vector< double > & | rStretches | |||
) | [inline, 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.
rDeformationGradients | A reference of a std::vector in which the deformation gradient in each element will be returned. Must be allocated prior to being passed in. | |
rStretches | A 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 399 of file AbstractCardiacMechanicsSolver.hpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mConstantFibreSheetDirections, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentElementFibreDirection, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mFibreSheetDirectionsDefinedByQuadraturePoint, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpCurrentElementFibreSheetMatrix, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpVariableFibreSheetDirections, and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::NUM_VERTICES_PER_ELEMENT.
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative | ( | AbstractMaterialLaw< DIM > * | pMaterialLaw, | |
c_matrix< double, DIM, DIM > & | rC, | |||
c_matrix< double, DIM, DIM > & | rInvC, | |||
double | pressure, | |||
unsigned | elementIndex, | |||
unsigned | currentQuadPointGlobalIndex, | |||
c_matrix< double, DIM, DIM > & | rT, | |||
FourthOrderTensor< DIM, DIM, DIM, DIM > & | rDTdE, | |||
bool | computeDTdE | |||
) | [inline, protected] |
Overloaded ComputeStressAndStressDerivative(), which computes the passive part of the stress as normal but also calls on the contraction model to get the active stress and adds it on.
pMaterialLaw | The material law for this element | |
rC | The Lagrangian deformation tensor (F^T F) | |
rInvC | The inverse of C. Should be computed by the user. | |
pressure | The current pressure | |
elementIndex | Index of the current element | |
currentQuadPointGlobalIndex | The index (assuming an outer loop over elements and an inner loop over quadrature points), of the current quadrature point. | |
rT | The stress will be returned in this parameter | |
rDTdE | the stress derivative will be returned in this parameter, assuming the final parameter is true | |
computeDTdE | A boolean flag saying whether the stress derivative is required or not. |
Definition at line 321 of file AbstractCardiacMechanicsSolver.hpp.
References AbstractMaterialLaw< DIM >::ComputeStressAndStressDerivative(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::IsImplicitSolver(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mConstantFibreSheetDirections, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentElementFibreDirection, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentTime, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mFibreSheetDirectionsDefinedByQuadraturePoint, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mNextTime, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpCurrentElementFibreSheetMatrix, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpVariableFibreSheetDirections, and AbstractMaterialLaw< DIM >::SetChangeOfBasisMatrix().
virtual void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs | ( | double | currentFibreStretch, | |
unsigned | currentQuadPointGlobalIndex, | |||
bool | assembleJacobian, | |||
double & | rActiveTension, | |||
double & | rDerivActiveTensionWrtLambda, | |||
double & | rDerivActiveTensionWrtDLambdaDt | |||
) | [protected, pure virtual] |
Pure method called in AbstractCardiacMechanicsSolver::ComputeStressAndStressDerivative(), 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.
currentFibreStretch | The stretch in the fibre direction | |
currentQuadPointGlobalIndex | quadrature point the integrand is currently being evaluated at in AssembleOnElement | |
assembleJacobian | A bool stating whether to assemble the Jacobian matrix. | |
rActiveTension | The returned active tension | |
rDerivActiveTensionWrtLambda | The returned dT_dLam, derivative of active tension wrt stretch. Only should be set in implicit solvers | |
rDerivActiveTensionWrtDLambdaDt | The 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 >.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative().
virtual GaussianQuadratureRule<DIM>* AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetQuadratureRule | ( | ) | [inline, virtual] |
Get the quadrature rule used in the elements.
Implements AbstractCardiacMechanicsSolverInterface< DIM >.
Definition at line 168 of file AbstractCardiacMechanicsSolver.hpp.
unsigned AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetTotalNumQuadPoints | ( | ) | [inline, virtual] |
Get the total number of quad points in the mesh. Pure, implemented in concrete solver
Implements AbstractCardiacMechanicsSolverInterface< DIM >.
Definition at line 162 of file AbstractCardiacMechanicsSolver.hpp.
virtual bool AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::IsImplicitSolver | ( | ) | [protected, pure virtual] |
Whether the solver is implicit or not (ie 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 >.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative().
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetCalciumAndVoltage | ( | std::vector< double > & | rCalciumConcentrations, | |
std::vector< double > & | rVoltages | |||
) | [inline, 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.
rCalciumConcentrations | Reference to a vector of intracellular calcium concentrations at each quadrature point | |
rVoltages | Reference to a vector of voltages at each quadrature point |
Implements AbstractCardiacMechanicsSolverInterface< DIM >.
Definition at line 302 of file AbstractCardiacMechanicsSolver.hpp.
References AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mContractionModelSystems, and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mTotalQuadPoints.
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetConstantFibreSheetDirections | ( | const c_matrix< double, DIM, DIM > & | rFibreSheetMatrix | ) | [inline, 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)
rFibreSheetMatrix | The 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 500 of file AbstractCardiacMechanicsSolver.hpp.
References EXCEPTION, and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mConstantFibreSheetDirections.
void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetVariableFibreSheetDirections | ( | std::string | orthoFile, | |
bool | definedPerQuadraturePoint | |||
) | [inline, 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.
orthoFile | the file containing the fibre/sheet directions | |
definedPerQuadraturePoint | whether 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 469 of file AbstractCardiacMechanicsSolver.hpp.
References RelativeTo::ChasteSourceRoot, EXCEPTION, FibreReader< DIM >::GetNextFibreSheetAndNormalMatrix(), FibreReader< DIM >::GetNumLinesOfData(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mFibreSheetDirectionsDefinedByQuadraturePoint, AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpVariableFibreSheetDirections, and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mTotalQuadPoints.
virtual void AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve | ( | double | time, | |
double | nextTime, | |||
double | odeTimestep | |||
) | [pure virtual] |
Solve for the deformation, integrating the contraction model ODEs.
time | the current time | |
nextTime | the next time | |
odeTimestep | the ODE timestep |
Implements AbstractCardiacMechanicsSolverInterface< DIM >.
Implemented in ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >, and ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, 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 69 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::AbstractCardiacMechanicsSolver(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative(), and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetConstantFibreSheetDirections().
std::vector<AbstractContractionModel*> AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mContractionModelSystems [protected] |
Vector of contraction model (pointers). One for each quadrature point. Note the indexing: the i-th entry corresponds to the i-th global quad point, when looping over elements and then quad points
Definition at line 49 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs(), ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::InitialiseContractionModels(), ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::InitialiseContractionModels(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetCalciumAndVoltage(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve(), ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve(), ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::~ExplicitCardiacMechanicsSolver(), and ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::~ImplicitCardiacMechanicsSolver().
c_vector<double,DIM> AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentElementFibreDirection [protected] |
The fibre direction for the current element being assembled on
Definition at line 86 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement(), and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative().
double AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mCurrentTime [protected] |
Current time
Definition at line 62 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve(), and ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve().
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 81 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative(), and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetVariableFibreSheetDirections().
double AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mNextTime [protected] |
Time to which the solver has been asked to solve to
Definition at line 64 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve(), and ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve().
double AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mOdeTimestep [protected] |
Time used to integrate the contraction model
Definition at line 66 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve(), and ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve().
c_matrix<double,DIM,DIM>* AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mpCurrentElementFibreSheetMatrix [protected] |
(Pointer to) the fibre-sheet matrix for the current element being assembled on
Definition at line 84 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement(), and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative().
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 75 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::AbstractCardiacMechanicsSolver(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeStressAndStressDerivative(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetVariableFibreSheetDirections(), and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::~AbstractCardiacMechanicsSolver().
std::vector<double> AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mStretches [protected] |
Stored stretches (in fibre direction, at each quadrature point). Should be stored when GetActiveTensionAndTensionDerivs() is called, and can be used either in that timestep (implicit solver), or the next timestep (explicit solver)
Definition at line 56 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::AbstractCardiacMechanicsSolver(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs(), ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::rGetFibreStretches(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve(), and ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::Solve().
unsigned AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mTotalQuadPoints [protected] |
Total number of quad points in the (mechanics) mesh
Definition at line 59 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::AbstractCardiacMechanicsSolver(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ImplicitCardiacMechanicsSolver(), ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::InitialiseContractionModels(), ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::InitialiseContractionModels(), AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetCalciumAndVoltage(), and AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::SetVariableFibreSheetDirections().
const unsigned AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::NUM_VERTICES_PER_ELEMENT = ELASTICITY_SOLVER::NUM_VERTICES_PER_ELEMENT [static, protected] |
Useful const from base class
Definition at line 42 of file AbstractCardiacMechanicsSolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::ComputeDeformationGradientAndStretchInEachElement().