ImplicitCardiacMechanicsSolver< DIM > Class Template Reference

#include <ImplicitCardiacMechanicsSolver.hpp>

Inheritance diagram for ImplicitCardiacMechanicsSolver< DIM >:

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

Collaboration graph
[legend]

List of all members.

Public Member Functions

 ImplicitCardiacMechanicsSolver (ContractionModel contractionModel, QuadraticMesh< DIM > *pQuadMesh, std::string outputDirectory, std::vector< unsigned > &rFixedNodes, AbstractIncompressibleMaterialLaw< DIM > *pMaterialLaw=NULL)
virtual ~ImplicitCardiacMechanicsSolver ()
std::vector< double > & rGetFibreStretches ()
void Solve (double time, double nextTime, double odeTimestep)

Private Member Functions

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

Private Attributes

std::vector< double > mStretchesLastTimeStep

Friends

class TestImplicitCardiacMechanicsSolver


Detailed Description

template<unsigned DIM>
class ImplicitCardiacMechanicsSolver< DIM >

Implicit Cardiac Mechanics Solver

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 51 of file ImplicitCardiacMechanicsSolver.hpp.


Constructor & Destructor Documentation

template<unsigned DIM>
ImplicitCardiacMechanicsSolver< DIM >::ImplicitCardiacMechanicsSolver ( ContractionModel  contractionModel,
QuadraticMesh< DIM > *  pQuadMesh,
std::string  outputDirectory,
std::vector< unsigned > &  rFixedNodes,
AbstractIncompressibleMaterialLaw< DIM > *  pMaterialLaw = NULL 
) [inline]

Constructor

Parameters:
contractionModel The contraction model.
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. Defaults to NULL, in which case a default material law is used.

Definition at line 35 of file ImplicitCardiacMechanicsSolver.cpp.

References EXCEPTION, AbstractCardiacMechanicsSolver< DIM >::mContractionModelSystems, ImplicitCardiacMechanicsSolver< DIM >::mStretchesLastTimeStep, and AbstractCardiacMechanicsSolver< DIM >::mTotalQuadPoints.

template<unsigned DIM>
ImplicitCardiacMechanicsSolver< DIM >::~ImplicitCardiacMechanicsSolver (  )  [inline, virtual]


Member Function Documentation

template<unsigned DIM>
bool ImplicitCardiacMechanicsSolver< DIM >::IsImplicitSolver (  )  [inline, private, virtual]

The current stretch ratio (in the fibre direction). Note the indexing: the i-th entry corresponds to the i-th global quad point, when looping over elements and then quad points This solver is an implicit solver (overloaded pure method)

Implements AbstractCardiacMechanicsSolver< DIM >.

Definition at line 69 of file ImplicitCardiacMechanicsSolver.hpp.

template<unsigned DIM>
void ImplicitCardiacMechanicsSolver< DIM >::GetActiveTensionAndTensionDerivs ( double  currentFibreStretch,
unsigned  currentQuadPointGlobalIndex,
bool  assembleJacobian,
double &  rActiveTension,
double &  rDerivActiveTensionWrtLambda,
double &  rDerivActiveTensionWrtDLambdaDt 
) [inline, private, virtual]

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:
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
rDerivActiveTensionWrtDLambdaDt The returned dT_dLamDot, derivative of active tension wrt stretch rate

Implements AbstractCardiacMechanicsSolver< DIM >.

Definition at line 136 of file ImplicitCardiacMechanicsSolver.cpp.

References EXCEPTION, AbstractContractionModel::GetNextActiveTension(), AbstractCardiacMechanicsSolver< DIM >::mContractionModelSystems, AbstractCardiacMechanicsSolver< DIM >::mCurrentTime, AbstractCardiacMechanicsSolver< DIM >::mNextTime, AbstractCardiacMechanicsSolver< DIM >::mOdeTimestep, AbstractCardiacMechanicsSolver< DIM >::mStretches, ImplicitCardiacMechanicsSolver< DIM >::mStretchesLastTimeStep, AbstractContractionModel::RunDoNotUpdate(), and AbstractContractionModel::SetStretchAndStretchRate().

template<unsigned DIM>
std::vector< double > & ImplicitCardiacMechanicsSolver< DIM >::rGetFibreStretches (  )  [inline]

Get lambda (the stretch ratio). NOTE: the i-th entry of this vector is assumed to be the i-th quad point obtained by looping over cells in the obvious way and then looping over quad points. These quad points, in the same order, can be obtained by using the QuadraturePointsGroup class.

Definition at line 98 of file ImplicitCardiacMechanicsSolver.cpp.

References AbstractCardiacMechanicsSolver< DIM >::mStretches.

template<unsigned DIM>
void ImplicitCardiacMechanicsSolver< DIM >::Solve ( double  time,
double  nextTime,
double  odeTimestep 
) [inline, 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:
time the current time
nextTime the next time
odeTimestep the ODE timestep

Implements AbstractCardiacMechanicsSolver< DIM >.

Definition at line 105 of file ImplicitCardiacMechanicsSolver.cpp.

References NonlinearElasticitySolver< DIM >::AssembleSystem(), AbstractNonlinearElasticitySolver< INCOMPRESSIBLE, DIM >::GetNumNewtonIterations(), AbstractCardiacMechanicsSolver< DIM >::mContractionModelSystems, AbstractCardiacMechanicsSolver< DIM >::mCurrentTime, AbstractCardiacMechanicsSolver< DIM >::mNextTime, AbstractCardiacMechanicsSolver< DIM >::mOdeTimestep, AbstractCardiacMechanicsSolver< DIM >::mStretches, ImplicitCardiacMechanicsSolver< DIM >::mStretchesLastTimeStep, and AbstractNonlinearElasticitySolver< INCOMPRESSIBLE, DIM >::Solve().


Member Data Documentation

template<unsigned DIM>
std::vector<double> ImplicitCardiacMechanicsSolver< DIM >::mStretchesLastTimeStep [private]

The stretch in the fibre direction at the last timestep, in order to compute the stretch rate. 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 61 of file ImplicitCardiacMechanicsSolver.hpp.

Referenced by ImplicitCardiacMechanicsSolver< DIM >::GetActiveTensionAndTensionDerivs(), ImplicitCardiacMechanicsSolver< DIM >::ImplicitCardiacMechanicsSolver(), and ImplicitCardiacMechanicsSolver< DIM >::Solve().


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

Generated on Mon Apr 18 11:37:17 2011 for Chaste by  doxygen 1.5.5