Chaste Commit::1fd4e48e3990e67db148bc1bc4cf6991a0049d0c
|
#include <AbstractCardiacMechanicsSolver.hpp>
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, DataAtQuadraturePoint > | mQuadPointToDataAtQuadPointMap |
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 |
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.
AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::AbstractCardiacMechanicsSolver | ( | QuadraticMesh< DIM > & | rQuadMesh, |
ElectroMechanicsProblemDefinition< DIM > & | rProblemDefinition, | ||
std::string | outputDirectory | ||
) |
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 41 of file AbstractCardiacMechanicsSolver.cpp.
|
virtual |
Destructor
Definition at line 145 of file AbstractCardiacMechanicsSolver.cpp.
|
protected |
Overloaded AddActiveStressAndStressDerivative(), which calls on the contraction model to get the active stress and add it on to the stress tensor
rC | The Lagrangian deformation tensor (F^T F) |
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 to be added to |
rDTdE | the stress derivative to be added to, assuming the final parameter is true |
addToDTdE | A boolean flag saying whether the stress derivative is required or not. |
Definition at line 208 of file AbstractCardiacMechanicsSolver.cpp.
References Determinant(), and Inverse().
|
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 418 of file AbstractCardiacMechanicsSolver.cpp.
References LinearBasisFunction< ELEMENT_DIM >::ComputeTransformedBasisFunctionDerivatives(), and AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex().
|
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.
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 >.
|
inlinevirtual |
Implements AbstractCardiacMechanicsSolverInterface< DIM >.
Definition at line 247 of file AbstractCardiacMechanicsSolver.hpp.
|
inlinevirtual |
Implements AbstractCardiacMechanicsSolverInterface< DIM >.
Definition at line 241 of file AbstractCardiacMechanicsSolver.hpp.
References AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mTotalQuadPoints.
|
protectedvirtual |
Sets relevant data at all quadrature points, including whether it is an active region or not.
Calls a contraction cell factory to assign a model to each (quadrature point in each) element.
Implements AbstractCardiacMechanicsSolverInterface< DIM >.
Definition at line 56 of file AbstractCardiacMechanicsSolver.cpp.
References DataAtQuadraturePoint_::ContractionModel, AbstractContractionCellFactory< DIM >::CreateContractionCellForElement(), EXCEPTION, AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetOwnership(), HeartRegionCode::GetValidBathId(), DataAtQuadraturePoint_::Stretch, and DataAtQuadraturePoint_::StretchLastTimeStep.
|
protectedpure virtual |
Implemented in ExplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >, and ImplicitCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >.
|
inline |
Definition at line 255 of file AbstractCardiacMechanicsSolver.hpp.
References AbstractCardiacMechanicsSolver< ELASTICITY_SOLVER, DIM >::mQuadPointToDataAtQuadPointMap.
|
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 165 of file AbstractCardiacMechanicsSolver.cpp.
References ContractionModelInputParameters_::intracellularCalciumConcentration, and ContractionModelInputParameters_::voltage.
|
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 517 of file AbstractCardiacMechanicsSolver.cpp.
References EXCEPTION.
|
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.
pMeshPair | the 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().
|
protected |
Over-ridden method which sets up an internal variable in the parent class, using the provided fibre-sheet direction information.
elementIndex | element global index |
currentQuadPointGlobalIndex | quad point global index |
Definition at line 190 of file AbstractCardiacMechanicsSolver.cpp.
|
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.
rOrthoFile | 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 487 of file AbstractCardiacMechanicsSolver.cpp.
References EXCEPTION, FileFinder::GetAbsolutePath(), FibreReader< DIM >::GetFibreSheetAndNormalMatrix(), and FibreReader< DIM >::GetNumLinesOfData().
|
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 >.
|
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.
|
protected |
The fibre direction for the current element being assembled on
Definition at line 137 of file AbstractCardiacMechanicsSolver.hpp.
|
protected |
The sheet direction for the current element being assembled on
Definition at line 140 of file AbstractCardiacMechanicsSolver.hpp.
|
protected |
The sheet normal direction for the current element being assembled on
Definition at line 143 of file AbstractCardiacMechanicsSolver.hpp.
|
protected |
Current time
Definition at line 110 of file AbstractCardiacMechanicsSolver.hpp.
|
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.
|
protected |
An iterator to the map, used to avoid having to repeatedly search the map.
Definition at line 101 of file AbstractCardiacMechanicsSolver.hpp.
|
protected |
Time to which the solver has been asked to solve to
Definition at line 113 of file AbstractCardiacMechanicsSolver.hpp.
|
protected |
Time used to integrate the contraction model
Definition at line 116 of file AbstractCardiacMechanicsSolver.hpp.
|
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.
|
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.
|
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().
|
protected |
This class contains all the information about the electro mechanics problem (except the material law)
Definition at line 149 of file AbstractCardiacMechanicsSolver.hpp.
|
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().
|
staticprotected |
Useful const from base class
Definition at line 84 of file AbstractCardiacMechanicsSolver.hpp.