AbstractCardiacMechanicsSolver< DIM > Class Template Reference

#include <AbstractCardiacMechanicsSolver.hpp>

Inheritance diagram for AbstractCardiacMechanicsSolver< DIM >:

Inheritance graph
[legend]
Collaboration diagram for AbstractCardiacMechanicsSolver< DIM >:

Collaboration graph
[legend]

List of all members.

Public Member Functions

 AbstractCardiacMechanicsSolver (QuadraticMesh< DIM > *pQuadMesh, std::string outputDirectory, std::vector< unsigned > &rFixedNodes, AbstractIncompressibleMaterialLaw< DIM > *pMaterialLaw)
 ~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 (AbstractIncompressibleMaterialLaw< 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
bool mAllocatedMaterialLawMemory
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 STENCIL_SIZE = NonlinearElasticitySolver<DIM>::STENCIL_SIZE
static const unsigned NUM_NODES_PER_ELEMENT = NonlinearElasticitySolver<DIM>::NUM_NODES_PER_ELEMENT
static const unsigned NUM_VERTICES_PER_ELEMENT = NonlinearElasticitySolver<DIM>::NUM_VERTICES_PER_ELEMENT


Detailed Description

template<unsigned DIM>
class AbstractCardiacMechanicsSolver< DIM >

AbstractCardiacMechanicsSolver

Base class to implicit and explicit cardiac mechanics solvers. Inherits from NonlinearElasticitySolver 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 49 of file AbstractCardiacMechanicsSolver.hpp.


Constructor & Destructor Documentation

template<unsigned DIM>
AbstractCardiacMechanicsSolver< DIM >::AbstractCardiacMechanicsSolver ( QuadraticMesh< DIM > *  pQuadMesh,
std::string  outputDirectory,
std::vector< unsigned > &  rFixedNodes,
AbstractIncompressibleMaterialLaw< DIM > *  pMaterialLaw 
) [inline]

template<unsigned DIM>
AbstractCardiacMechanicsSolver< DIM >::~AbstractCardiacMechanicsSolver (  )  [inline]


Member Function Documentation

template<unsigned DIM>
virtual bool AbstractCardiacMechanicsSolver< 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< DIM >, and ImplicitCardiacMechanicsSolver< DIM >.

Referenced by AbstractCardiacMechanicsSolver< DIM >::ComputeStressAndStressDerivative().

template<unsigned DIM>
void AbstractCardiacMechanicsSolver< DIM >::ComputeStressAndStressDerivative ( AbstractIncompressibleMaterialLaw< 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, virtual]

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.

Parameters:
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.

Reimplemented from NonlinearElasticitySolver< DIM >.

Definition at line 339 of file AbstractCardiacMechanicsSolver.hpp.

References AbstractIncompressibleMaterialLaw< DIM >::ComputeStressAndStressDerivative(), AbstractCardiacMechanicsSolver< DIM >::GetActiveTensionAndTensionDerivs(), AbstractCardiacMechanicsSolver< DIM >::IsImplicitSolver(), AbstractCardiacMechanicsSolver< DIM >::mConstantFibreSheetDirections, AbstractCardiacMechanicsSolver< DIM >::mCurrentElementFibreDirection, AbstractCardiacMechanicsSolver< DIM >::mCurrentTime, AbstractCardiacMechanicsSolver< DIM >::mFibreSheetDirectionsDefinedByQuadraturePoint, AbstractCardiacMechanicsSolver< DIM >::mNextTime, AbstractCardiacMechanicsSolver< DIM >::mpCurrentElementFibreSheetMatrix, AbstractCardiacMechanicsSolver< DIM >::mpVariableFibreSheetDirections, and AbstractIncompressibleMaterialLaw< DIM >::SetChangeOfBasisMatrix().

template<unsigned DIM>
virtual void AbstractCardiacMechanicsSolver< 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.

Parameters:
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< DIM >, and ImplicitCardiacMechanicsSolver< DIM >.

Referenced by AbstractCardiacMechanicsSolver< DIM >::ComputeStressAndStressDerivative().

template<unsigned DIM>
unsigned AbstractCardiacMechanicsSolver< DIM >::GetTotalNumQuadPoints (  )  [inline]

Get the total number of quad points in the mesh. Pure, implemented in concrete solver

Definition at line 181 of file AbstractCardiacMechanicsSolver.hpp.

References AbstractCardiacMechanicsSolver< DIM >::mTotalQuadPoints.

template<unsigned DIM>
virtual GaussianQuadratureRule<DIM>* AbstractCardiacMechanicsSolver< DIM >::GetQuadratureRule (  )  [inline, virtual]

Get the quadrature rule used in the elements.

Definition at line 187 of file AbstractCardiacMechanicsSolver.hpp.

References NonlinearElasticitySolver< DIM >::mpQuadratureRule.

template<unsigned DIM>
void AbstractCardiacMechanicsSolver< 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)

Parameters:
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 522 of file AbstractCardiacMechanicsSolver.hpp.

References EXCEPTION, and AbstractCardiacMechanicsSolver< DIM >::mConstantFibreSheetDirections.

template<unsigned DIM>
void AbstractCardiacMechanicsSolver< DIM >::SetVariableFibreSheetDirections ( std::string  orthoFile,
bool  definedPerQuadraturePoint 
) [inline]

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:
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).

Definition at line 487 of file AbstractCardiacMechanicsSolver.hpp.

References RelativeTo::ChasteSourceRoot, EXCEPTION, FibreReader< DIM >::GetNextFibreSheetAndNormalMatrix(), FibreReader< DIM >::GetNumLinesOfData(), AbstractCardiacMechanicsSolver< DIM >::mFibreSheetDirectionsDefinedByQuadraturePoint, NonlinearElasticitySolver< DIM >::mpQuadMesh, AbstractCardiacMechanicsSolver< DIM >::mpVariableFibreSheetDirections, and AbstractCardiacMechanicsSolver< DIM >::mTotalQuadPoints.

template<unsigned DIM>
void AbstractCardiacMechanicsSolver< 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.

Parameters:
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 320 of file AbstractCardiacMechanicsSolver.hpp.

References AbstractCardiacMechanicsSolver< DIM >::mContractionModelSystems, and AbstractCardiacMechanicsSolver< DIM >::mTotalQuadPoints.

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

Solve for the deformation, integrating the contraction model ODEs.

Parameters:
time the current time
nextTime the next time
odeTimestep the ODE timestep

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

template<unsigned DIM>
void AbstractCardiacMechanicsSolver< DIM >::ComputeDeformationGradientAndStretchInEachElement ( std::vector< c_matrix< double, DIM, DIM > > &  rDeformationGradients,
std::vector< double > &  rStretches 
) [inline]

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:
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.

Definition at line 417 of file AbstractCardiacMechanicsSolver.hpp.

References LinearBasisFunction< ELEMENT_DIM >::ComputeTransformedBasisFunctionDerivatives(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractCardiacMechanicsSolver< DIM >::mConstantFibreSheetDirections, AbstractCardiacMechanicsSolver< DIM >::mCurrentElementFibreDirection, AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, AbstractCardiacMechanicsSolver< DIM >::mFibreSheetDirectionsDefinedByQuadraturePoint, AbstractCardiacMechanicsSolver< DIM >::mpCurrentElementFibreSheetMatrix, NonlinearElasticitySolver< DIM >::mpQuadMesh, AbstractCardiacMechanicsSolver< DIM >::mpVariableFibreSheetDirections, and AbstractCardiacMechanicsSolver< DIM >::NUM_VERTICES_PER_ELEMENT.


Member Data Documentation

template<unsigned DIM>
const unsigned AbstractCardiacMechanicsSolver< DIM >::STENCIL_SIZE = NonlinearElasticitySolver<DIM>::STENCIL_SIZE [static, protected]

Stencil size

Reimplemented from NonlinearElasticitySolver< DIM >.

Definition at line 52 of file AbstractCardiacMechanicsSolver.hpp.

template<unsigned DIM>
const unsigned AbstractCardiacMechanicsSolver< DIM >::NUM_NODES_PER_ELEMENT = NonlinearElasticitySolver<DIM>::NUM_NODES_PER_ELEMENT [static, protected]

Number of nodes per element

Reimplemented from NonlinearElasticitySolver< DIM >.

Definition at line 53 of file AbstractCardiacMechanicsSolver.hpp.

template<unsigned DIM>
const unsigned AbstractCardiacMechanicsSolver< DIM >::NUM_VERTICES_PER_ELEMENT = NonlinearElasticitySolver<DIM>::NUM_VERTICES_PER_ELEMENT [static, protected]

template<unsigned DIM>
std::vector<AbstractContractionModel*> AbstractCardiacMechanicsSolver< DIM >::mContractionModelSystems [protected]

template<unsigned DIM>
std::vector<double> AbstractCardiacMechanicsSolver< DIM >::mStretches [protected]

template<unsigned DIM>
unsigned AbstractCardiacMechanicsSolver< DIM >::mTotalQuadPoints [protected]

template<unsigned DIM>
bool AbstractCardiacMechanicsSolver< DIM >::mAllocatedMaterialLawMemory [protected]

template<unsigned DIM>
double AbstractCardiacMechanicsSolver< DIM >::mCurrentTime [protected]

template<unsigned DIM>
double AbstractCardiacMechanicsSolver< DIM >::mNextTime [protected]

template<unsigned DIM>
double AbstractCardiacMechanicsSolver< DIM >::mOdeTimestep [protected]

template<unsigned DIM>
c_matrix<double,DIM,DIM> AbstractCardiacMechanicsSolver< DIM >::mConstantFibreSheetDirections [protected]

template<unsigned DIM>
std::vector<c_matrix<double,DIM,DIM> >* AbstractCardiacMechanicsSolver< DIM >::mpVariableFibreSheetDirections [protected]

template<unsigned DIM>
bool AbstractCardiacMechanicsSolver< 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 96 of file AbstractCardiacMechanicsSolver.hpp.

Referenced by AbstractCardiacMechanicsSolver< DIM >::ComputeDeformationGradientAndStretchInEachElement(), AbstractCardiacMechanicsSolver< DIM >::ComputeStressAndStressDerivative(), and AbstractCardiacMechanicsSolver< DIM >::SetVariableFibreSheetDirections().

template<unsigned DIM>
c_matrix<double,DIM,DIM>* AbstractCardiacMechanicsSolver< DIM >::mpCurrentElementFibreSheetMatrix [protected]

template<unsigned DIM>
c_vector<double,DIM> AbstractCardiacMechanicsSolver< DIM >::mCurrentElementFibreDirection [protected]


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

Generated on Mon Nov 1 12:35:33 2010 for Chaste by  doxygen 1.5.5