#include <AbstractCardiacMechanicsAssembler.hpp>
Inherits NonlinearElasticityAssembler< DIM >.
Inherited by ExplicitCardiacMechanicsAssembler< DIM >, and ImplicitCardiacMechanicsAssembler< DIM >.
Public Member Functions | |
AbstractCardiacMechanicsAssembler (QuadraticMesh< DIM > *pQuadMesh, std::string outputDirectory, std::vector< unsigned > &rFixedNodes, AbstractIncompressibleMaterialLaw< DIM > *pMaterialLaw) | |
~AbstractCardiacMechanicsAssembler () | |
unsigned | GetTotalNumQuadPoints () |
virtual GaussianQuadratureRule < DIM > * | GetQuadratureRule () |
void | SetConstantFibreSheetDirections (const c_matrix< double, DIM, DIM > &rFibreSheetMatrix) |
void | SetVariableFibreSheetDirections (std::string orthoFile) |
void | SetCalciumAndVoltage (std::vector< double > &rCalciumConcentrations, std::vector< double > &rVoltages) |
virtual void | Solve (double time, double nextTime, double odeTimestep)=0 |
Protected Member Functions | |
void | CheckOrthogonality (c_matrix< double, DIM, DIM > &rMatrix) |
virtual bool | IsImplicitSolver ()=0 |
void | AssembleOnElement (Element< DIM, DIM > &rElement, c_matrix< double, STENCIL_SIZE, STENCIL_SIZE > &rAElem, c_matrix< double, STENCIL_SIZE, STENCIL_SIZE > &rAElemPrecond, c_vector< double, STENCIL_SIZE > &rBElem, bool assembleResidual, bool assembleJacobian) |
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 |
bool | mAllocatedMaterialLawMemory |
double | mCurrentTime |
double | mNextTime |
double | mOdeTimestep |
c_matrix< double, DIM, DIM > | mConstantFibreSheetDirections |
std::vector< c_matrix< double, DIM, DIM > > * | mpVariableFibreSheetDirections |
Static Protected Attributes | |
static const unsigned | STENCIL_SIZE = NonlinearElasticityAssembler<DIM>::STENCIL_SIZE |
static const unsigned | NUM_NODES_PER_ELEMENT = NonlinearElasticityAssembler<DIM>::NUM_NODES_PER_ELEMENT |
static const unsigned | NUM_VERTICES_PER_ELEMENT = NonlinearElasticityAssembler<DIM>::NUM_VERTICES_PER_ELEMENT |
AbstractCardiacMechanicsAssembler
Base class to implicit and explicit cardiac mechanics assemblers. Inherits from NonlinearElasticityAssembler Main method is the overloaded AssembleOnElement which does the extra work needed for cardiac problems. The child classes hold the contraction models and need to implement a method for getting the active tension from the model.
Definition at line 47 of file AbstractCardiacMechanicsAssembler.hpp.
AbstractCardiacMechanicsAssembler< DIM >::AbstractCardiacMechanicsAssembler | ( | QuadraticMesh< DIM > * | pQuadMesh, | |
std::string | outputDirectory, | |||
std::vector< unsigned > & | rFixedNodes, | |||
AbstractIncompressibleMaterialLaw< DIM > * | pMaterialLaw | |||
) | [inline] |
Constructor
pQuadMesh | A pointer to the mesh. | |
outputDirectory | The output directory, relative to TEST_OUTPUT | |
rFixedNodes | The fixed nodes | |
pMaterialLaw | The material law for the tissue. If NULL the default (NashHunterPoleZero) law is used. |
Definition at line 177 of file AbstractCardiacMechanicsAssembler.hpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), GaussianQuadratureRule< ELEMENT_DIM >::GetNumQuadPoints(), AbstractCardiacMechanicsAssembler< DIM >::mAllocatedMaterialLawMemory, AbstractCardiacMechanicsAssembler< DIM >::mConstantFibreSheetDirections, NonlinearElasticityAssembler< DIM >::mpQuadratureRule, AbstractCardiacMechanicsAssembler< DIM >::mpVariableFibreSheetDirections, AbstractCardiacMechanicsAssembler< DIM >::mStretches, and AbstractCardiacMechanicsAssembler< DIM >::mTotalQuadPoints.
AbstractCardiacMechanicsAssembler< DIM >::~AbstractCardiacMechanicsAssembler | ( | ) | [inline] |
Destructor just deletes memory if it was allocated
Definition at line 214 of file AbstractCardiacMechanicsAssembler.hpp.
References AbstractCardiacMechanicsAssembler< DIM >::mAllocatedMaterialLawMemory, AbstractNonlinearElasticityAssembler< DIM >::mMaterialLaws, and AbstractCardiacMechanicsAssembler< DIM >::mpVariableFibreSheetDirections.
void AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement | ( | Element< DIM, DIM > & | rElement, | |
c_matrix< double, STENCIL_SIZE, STENCIL_SIZE > & | rAElem, | |||
c_matrix< double, STENCIL_SIZE, STENCIL_SIZE > & | rAElemPrecond, | |||
c_vector< double, STENCIL_SIZE > & | rBElem, | |||
bool | assembleResidual, | |||
bool | assembleJacobian | |||
) | [inline, protected, virtual] |
Overloaded AssembleOnElement. Apart from a tiny bit of initial set up and the lack of the body force term in the residual, the bits where this is different to the base class AssembleOnElement are restricted to two bits (see code): getting Ta and using it to compute the stress, and (in when Ta is a function of the stretch) the addition of extra term to the Jacobian.
rElement | The element to assemble on. | |
rAElem | The element's contribution to the LHS matrix is returned in this n by n matrix, where n is the no. of nodes in this element. There is no need to zero this matrix before calling. | |
rAElemPrecond | The element's contribution to the matrix passed to PetSC in creating a preconditioner | |
rBElem | The element's contribution to the RHS vector is returned in this vector of length n, the no. of nodes in this element. There is no need to zero this vector before calling. | |
assembleResidual | A bool stating whether to assemble the residual vector. | |
assembleJacobian | A bool stating whether to assemble the Jacobian matrix. |
Reimplemented from NonlinearElasticityAssembler< DIM >.
Definition at line 312 of file AbstractCardiacMechanicsAssembler.hpp.
References AbstractIncompressibleMaterialLaw< DIM >::ComputeStressAndStressDerivative(), AbstractNonlinearElasticityAssembler< DIM >::dTdE, AbstractCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), GaussianQuadratureRule< ELEMENT_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEMENT_DIM >::GetWeight(), AbstractCardiacMechanicsAssembler< DIM >::IsImplicitSolver(), AbstractCardiacMechanicsAssembler< DIM >::mConstantFibreSheetDirections, AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, AbstractCardiacMechanicsAssembler< DIM >::mCurrentTime, AbstractNonlinearElasticityAssembler< DIM >::mMaterialLaws, AbstractCardiacMechanicsAssembler< DIM >::mNextTime, AbstractCardiacMechanicsAssembler< DIM >::mOdeTimestep, NonlinearElasticityAssembler< DIM >::mpQuadMesh, NonlinearElasticityAssembler< DIM >::mpQuadratureRule, AbstractCardiacMechanicsAssembler< DIM >::mpVariableFibreSheetDirections, AbstractCardiacMechanicsAssembler< DIM >::NUM_NODES_PER_ELEMENT, AbstractCardiacMechanicsAssembler< DIM >::NUM_VERTICES_PER_ELEMENT, GaussianQuadratureRule< ELEMENT_DIM >::rGetQuadPoint(), FourthOrderTensor< DIM >::SetAsProduct(), AbstractIncompressibleMaterialLaw< DIM >::SetChangeOfBasisMatrix(), and AbstractCardiacMechanicsAssembler< DIM >::STENCIL_SIZE.
void AbstractCardiacMechanicsAssembler< DIM >::CheckOrthogonality | ( | c_matrix< double, DIM, DIM > & | rMatrix | ) | [inline, protected] |
Check whether the given matrix is orthogonal, by computing A^T A and verifying whether each component matches the identity matrix, to a given tolerance.
rMatrix | reference to the matrix being tested |
Definition at line 94 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::SetConstantFibreSheetDirections(), and AbstractCardiacMechanicsAssembler< DIM >::SetVariableFibreSheetDirections().
virtual void AbstractCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs | ( | double | currentFibreStretch, | |
unsigned | currentQuadPointGlobalIndex, | |||
bool | assembleJacobian, | |||
double & | rActiveTension, | |||
double & | rDerivActiveTensionWrtLambda, | |||
double & | rDerivActiveTensionWrtDLambdaDt | |||
) | [protected, pure virtual] |
Pure method called in AbstractCardiacMechanicsAssembler::AssembleOnElement(), 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 integrand 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 assemblers | |
rDerivActiveTensionWrtDLambdaDt | The returned dT_dLamDot, derivative of active tension wrt stretch rate. Only should be set in implicit assemblers |
Implemented in ExplicitCardiacMechanicsAssembler< DIM >, and ImplicitCardiacMechanicsAssembler< DIM >.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement().
virtual GaussianQuadratureRule<DIM>* AbstractCardiacMechanicsAssembler< DIM >::GetQuadratureRule | ( | ) | [inline, virtual] |
Get the quadrature rule used in the elements.
Definition at line 236 of file AbstractCardiacMechanicsAssembler.hpp.
References NonlinearElasticityAssembler< DIM >::mpQuadratureRule.
unsigned AbstractCardiacMechanicsAssembler< DIM >::GetTotalNumQuadPoints | ( | ) | [inline] |
Get the total number of quad points in the mesh. Pure, implemented in concrete assembler
Definition at line 230 of file AbstractCardiacMechanicsAssembler.hpp.
References AbstractCardiacMechanicsAssembler< DIM >::mTotalQuadPoints.
virtual bool AbstractCardiacMechanicsAssembler< 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 ExplicitCardiacMechanicsAssembler< DIM >, and ImplicitCardiacMechanicsAssembler< DIM >.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement().
void AbstractCardiacMechanicsAssembler< DIM >::SetCalciumAndVoltage | ( | std::vector< double > & | rCalciumConcentrations, | |
std::vector< double > & | rVoltages | |||
) | [inline] |
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 |
Definition at line 292 of file AbstractCardiacMechanicsAssembler.hpp.
References AbstractCardiacMechanicsAssembler< DIM >::mContractionModelSystems, and AbstractCardiacMechanicsAssembler< DIM >::mTotalQuadPoints.
void AbstractCardiacMechanicsAssembler< DIM >::SetConstantFibreSheetDirections | ( | const c_matrix< double, DIM, DIM > & | rFibreSheetMatrix | ) | [inline] |
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). |
Definition at line 248 of file AbstractCardiacMechanicsAssembler.hpp.
References AbstractCardiacMechanicsAssembler< DIM >::CheckOrthogonality(), and AbstractCardiacMechanicsAssembler< DIM >::mConstantFibreSheetDirections.
void AbstractCardiacMechanicsAssembler< DIM >::SetVariableFibreSheetDirections | ( | std::string | orthoFile | ) | [inline] |
Set a variable fibre-sheet-normal direction (matrices), one for each element, from a file. The file should be a .ortho file (ie each line has the fibre dir, sheet dir, normal dir for that element). The number of elements must match the number in the MECHANICS mesh!
orthoFile | the file containing the fibre/sheet directions |
Definition at line 680 of file AbstractCardiacMechanicsAssembler.hpp.
References AbstractCardiacMechanicsAssembler< DIM >::CheckOrthogonality(), NonlinearElasticityAssembler< DIM >::mpQuadMesh, and AbstractCardiacMechanicsAssembler< DIM >::mpVariableFibreSheetDirections.
virtual void AbstractCardiacMechanicsAssembler< 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 |
Implemented in ExplicitCardiacMechanicsAssembler< DIM >, and ImplicitCardiacMechanicsAssembler< DIM >.
bool AbstractCardiacMechanicsAssembler< DIM >::mAllocatedMaterialLawMemory [protected] |
Whether the material law was passed in or the default used
Definition at line 72 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AbstractCardiacMechanicsAssembler(), and AbstractCardiacMechanicsAssembler< DIM >::~AbstractCardiacMechanicsAssembler().
c_matrix<double,DIM,DIM> AbstractCardiacMechanicsAssembler< 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 82 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AbstractCardiacMechanicsAssembler(), AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), and AbstractCardiacMechanicsAssembler< DIM >::SetConstantFibreSheetDirections().
std::vector<AbstractContractionModel*> AbstractCardiacMechanicsAssembler< 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 59 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by ExplicitCardiacMechanicsAssembler< DIM >::ExplicitCardiacMechanicsAssembler(), ImplicitCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs(), ExplicitCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsAssembler< DIM >::ImplicitCardiacMechanicsAssembler(), AbstractCardiacMechanicsAssembler< DIM >::SetCalciumAndVoltage(), ImplicitCardiacMechanicsAssembler< DIM >::Solve(), ExplicitCardiacMechanicsAssembler< DIM >::Solve(), ExplicitCardiacMechanicsAssembler< DIM >::~ExplicitCardiacMechanicsAssembler(), and ImplicitCardiacMechanicsAssembler< DIM >::~ImplicitCardiacMechanicsAssembler().
double AbstractCardiacMechanicsAssembler< DIM >::mCurrentTime [protected] |
Current time
Definition at line 75 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), ImplicitCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsAssembler< DIM >::Solve(), and ExplicitCardiacMechanicsAssembler< DIM >::Solve().
double AbstractCardiacMechanicsAssembler< DIM >::mNextTime [protected] |
Time to which the solver has been asked to solve to
Definition at line 77 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), ImplicitCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsAssembler< DIM >::Solve(), and ExplicitCardiacMechanicsAssembler< DIM >::Solve().
double AbstractCardiacMechanicsAssembler< DIM >::mOdeTimestep [protected] |
Time used to integrate the contraction model
Definition at line 79 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), ImplicitCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsAssembler< DIM >::Solve(), and ExplicitCardiacMechanicsAssembler< DIM >::Solve().
std::vector<c_matrix<double,DIM,DIM> >* AbstractCardiacMechanicsAssembler< 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 88 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AbstractCardiacMechanicsAssembler(), AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), AbstractCardiacMechanicsAssembler< DIM >::SetVariableFibreSheetDirections(), and AbstractCardiacMechanicsAssembler< DIM >::~AbstractCardiacMechanicsAssembler().
std::vector<double> AbstractCardiacMechanicsAssembler< 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 66 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AbstractCardiacMechanicsAssembler(), ImplicitCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs(), ExplicitCardiacMechanicsAssembler< DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsAssembler< DIM >::rGetFibreStretches(), ImplicitCardiacMechanicsAssembler< DIM >::Solve(), and ExplicitCardiacMechanicsAssembler< DIM >::Solve().
unsigned AbstractCardiacMechanicsAssembler< DIM >::mTotalQuadPoints [protected] |
Total number of quad points in the (mechanics) mesh
Definition at line 69 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AbstractCardiacMechanicsAssembler(), ExplicitCardiacMechanicsAssembler< DIM >::ExplicitCardiacMechanicsAssembler(), AbstractCardiacMechanicsAssembler< DIM >::GetTotalNumQuadPoints(), ImplicitCardiacMechanicsAssembler< DIM >::ImplicitCardiacMechanicsAssembler(), and AbstractCardiacMechanicsAssembler< DIM >::SetCalciumAndVoltage().
const unsigned AbstractCardiacMechanicsAssembler< DIM >::NUM_NODES_PER_ELEMENT = NonlinearElasticityAssembler<DIM>::NUM_NODES_PER_ELEMENT [static, protected] |
Number of nodes per element
Reimplemented from NonlinearElasticityAssembler< DIM >.
Definition at line 51 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement().
const unsigned AbstractCardiacMechanicsAssembler< DIM >::NUM_VERTICES_PER_ELEMENT = NonlinearElasticityAssembler<DIM>::NUM_VERTICES_PER_ELEMENT [static, protected] |
Number of vertices per element
Reimplemented from NonlinearElasticityAssembler< DIM >.
Definition at line 52 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement().
const unsigned AbstractCardiacMechanicsAssembler< DIM >::STENCIL_SIZE = NonlinearElasticityAssembler<DIM>::STENCIL_SIZE [static, protected] |
Stencil size
Reimplemented from NonlinearElasticityAssembler< DIM >.
Definition at line 50 of file AbstractCardiacMechanicsAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement().