ImplicitCardiacMechanicsAssembler< DIM > Class Template Reference

#include <ImplicitCardiacMechanicsAssembler.hpp>

Inheritance diagram for ImplicitCardiacMechanicsAssembler< DIM >:

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

Collaboration graph
[legend]

List of all members.

Public Member Functions

 ImplicitCardiacMechanicsAssembler (QuadraticMesh< DIM > *pQuadMesh, std::string outputDirectory, std::vector< unsigned > &rFixedNodes, AbstractIncompressibleMaterialLaw< DIM > *pMaterialLaw=NULL)
 ~ImplicitCardiacMechanicsAssembler ()
unsigned GetTotalNumQuadPoints ()
GaussianQuadratureRule< DIM > * GetQuadratureRule ()
void SetIntracellularCalciumConcentrations (std::vector< double > &caI)
std::vector< double > & rGetLambda ()
void Solve (double time, double nextTime, double odeTimestep)

Private Member Functions

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)

Private Attributes

std::vector
< NhsSystemWithImplicitSolver
mCellMechSystems
std::vector< double > mLambdaLastTimeStep
std::vector< double > mLambda
double mCurrentTime
double mNextTime
double mOdeTimestep
bool mAllocatedMaterialLawMemory
unsigned mTotalQuadPoints

Static Private 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

Friends

class TestImplicitCardiacMechanicsAssembler


Detailed Description

template<unsigned DIM>
class ImplicitCardiacMechanicsAssembler< DIM >

Implicit Cardiac Mechanics Assembler

Solves cardiac mechanics implicitly (together with the NHS cell models for determining the active tension), taking in the intracellular Calcium concentration. See CardiacElectroMechanicsProblem documentation for more detail.

Definition at line 51 of file ImplicitCardiacMechanicsAssembler.hpp.


Constructor & Destructor Documentation

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

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


Member Function Documentation

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

Get the total number of quad points in the mesh

Definition at line 71 of file ImplicitCardiacMechanicsAssembler.cpp.

References ImplicitCardiacMechanicsAssembler< DIM >::mTotalQuadPoints.

template<unsigned DIM>
GaussianQuadratureRule< DIM > * ImplicitCardiacMechanicsAssembler< DIM >::GetQuadratureRule (  )  [inline]

Get the quadrature rule used in the elements

Definition at line 77 of file ImplicitCardiacMechanicsAssembler.cpp.

References NonlinearElasticityAssembler< DIM >::mpQuadratureRule.

template<unsigned DIM>
void ImplicitCardiacMechanicsAssembler< DIM >::SetIntracellularCalciumConcentrations ( std::vector< double > &  caI  )  [inline]

Set the intracellular Calcium concentrations (note: in an explicit algorithm we would set the active tension as the forcing quantity; the implicit algorithm takes in the Calcium concentration and solves for the active tension implicitly together with the mechanics).

Parameters:
caI the intracellular calcium concentrations

Definition at line 83 of file ImplicitCardiacMechanicsAssembler.cpp.

References ImplicitCardiacMechanicsAssembler< DIM >::mCellMechSystems, and ImplicitCardiacMechanicsAssembler< DIM >::mTotalQuadPoints.

template<unsigned DIM>
std::vector< double > & ImplicitCardiacMechanicsAssembler< DIM >::rGetLambda (  )  [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 93 of file ImplicitCardiacMechanicsAssembler.cpp.

References ImplicitCardiacMechanicsAssembler< DIM >::mLambda.

template<unsigned DIM>
void ImplicitCardiacMechanicsAssembler< DIM >::Solve ( double  time,
double  nextTime,
double  odeTimestep 
) [inline]

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

Definition at line 100 of file ImplicitCardiacMechanicsAssembler.cpp.

References NonlinearElasticityAssembler< DIM >::AssembleSystem(), ImplicitCardiacMechanicsAssembler< DIM >::mCellMechSystems, ImplicitCardiacMechanicsAssembler< DIM >::mCurrentTime, ImplicitCardiacMechanicsAssembler< DIM >::mLambdaLastTimeStep, ImplicitCardiacMechanicsAssembler< DIM >::mNextTime, ImplicitCardiacMechanicsAssembler< DIM >::mOdeTimestep, and AbstractNonlinearElasticityAssembler< DIM >::Solve().

template<unsigned DIM>
void ImplicitCardiacMechanicsAssembler< 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, private, 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): calculating Ta implicitly and using it to compute the stress, and the addition of a corresponding extra term to the Jacobian.

Parameters:
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 128 of file ImplicitCardiacMechanicsAssembler.cpp.

References QuadraticBasisFunction< ELEM_DIM >::ComputeBasisFunctions(), LinearBasisFunction< ELEM_DIM >::ComputeBasisFunctions(), AbstractIncompressibleMaterialLaw< DIM >::ComputeStressAndStressDerivative(), QuadraticBasisFunction< ELEM_DIM >::ComputeTransformedBasisFunctionDerivatives(), AbstractNonlinearElasticityAssembler< DIM >::dTdE, NhsSystemWithImplicitSolver::GetActiveTensionAtNextTime(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), GaussianQuadratureRule< ELEM_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEM_DIM >::GetWeight(), ImplicitCardiacMechanicsAssembler< DIM >::mCellMechSystems, AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, ImplicitCardiacMechanicsAssembler< DIM >::mCurrentTime, ImplicitCardiacMechanicsAssembler< DIM >::mLambda, ImplicitCardiacMechanicsAssembler< DIM >::mLambdaLastTimeStep, AbstractNonlinearElasticityAssembler< DIM >::mMaterialLaws, ImplicitCardiacMechanicsAssembler< DIM >::mNextTime, ImplicitCardiacMechanicsAssembler< DIM >::mOdeTimestep, NonlinearElasticityAssembler< DIM >::mpQuadMesh, NonlinearElasticityAssembler< DIM >::mpQuadratureRule, ImplicitCardiacMechanicsAssembler< DIM >::NUM_NODES_PER_ELEMENT, ImplicitCardiacMechanicsAssembler< DIM >::NUM_VERTICES_PER_ELEMENT, GaussianQuadratureRule< ELEM_DIM >::rGetQuadPoint(), NhsSystemWithImplicitSolver::SetActiveTensionInitialGuess(), FourthOrderTensor< DIM >::SetAsProduct(), NhsCellularMechanicsOdeSystem::SetLambdaAndDerivative(), and NhsSystemWithImplicitSolver::SolveDoNotUpdate().


Member Data Documentation

template<unsigned DIM>
const unsigned ImplicitCardiacMechanicsAssembler< DIM >::STENCIL_SIZE = NonlinearElasticityAssembler<DIM>::STENCIL_SIZE [static, private]

Stencil size

Reimplemented from NonlinearElasticityAssembler< DIM >.

Definition at line 56 of file ImplicitCardiacMechanicsAssembler.hpp.

template<unsigned DIM>
const unsigned ImplicitCardiacMechanicsAssembler< DIM >::NUM_NODES_PER_ELEMENT = NonlinearElasticityAssembler<DIM>::NUM_NODES_PER_ELEMENT [static, private]

Number of nodes per element

Reimplemented from NonlinearElasticityAssembler< DIM >.

Definition at line 57 of file ImplicitCardiacMechanicsAssembler.hpp.

Referenced by ImplicitCardiacMechanicsAssembler< DIM >::AssembleOnElement().

template<unsigned DIM>
const unsigned ImplicitCardiacMechanicsAssembler< DIM >::NUM_VERTICES_PER_ELEMENT = NonlinearElasticityAssembler<DIM>::NUM_VERTICES_PER_ELEMENT [static, private]

Number of vertices per element

Reimplemented from NonlinearElasticityAssembler< DIM >.

Definition at line 58 of file ImplicitCardiacMechanicsAssembler.hpp.

Referenced by ImplicitCardiacMechanicsAssembler< DIM >::AssembleOnElement().

template<unsigned DIM>
std::vector<NhsSystemWithImplicitSolver> ImplicitCardiacMechanicsAssembler< DIM >::mCellMechSystems [private]

The NHS cell systems (with their own implicit solvers, which take in [Ca]_i and return Ta. 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 64 of file ImplicitCardiacMechanicsAssembler.hpp.

Referenced by ImplicitCardiacMechanicsAssembler< DIM >::AssembleOnElement(), ImplicitCardiacMechanicsAssembler< DIM >::ImplicitCardiacMechanicsAssembler(), ImplicitCardiacMechanicsAssembler< DIM >::SetIntracellularCalciumConcentrations(), and ImplicitCardiacMechanicsAssembler< DIM >::Solve().

template<unsigned DIM>
std::vector<double> ImplicitCardiacMechanicsAssembler< DIM >::mLambdaLastTimeStep [private]

The stretch ratio (in the fibre direction) at the last timestep. 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 70 of file ImplicitCardiacMechanicsAssembler.hpp.

Referenced by ImplicitCardiacMechanicsAssembler< DIM >::AssembleOnElement(), ImplicitCardiacMechanicsAssembler< DIM >::ImplicitCardiacMechanicsAssembler(), and ImplicitCardiacMechanicsAssembler< DIM >::Solve().

template<unsigned DIM>
std::vector<double> ImplicitCardiacMechanicsAssembler< DIM >::mLambda [private]

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

Definition at line 76 of file ImplicitCardiacMechanicsAssembler.hpp.

Referenced by ImplicitCardiacMechanicsAssembler< DIM >::AssembleOnElement(), ImplicitCardiacMechanicsAssembler< DIM >::ImplicitCardiacMechanicsAssembler(), and ImplicitCardiacMechanicsAssembler< DIM >::rGetLambda().

template<unsigned DIM>
double ImplicitCardiacMechanicsAssembler< DIM >::mCurrentTime [private]

template<unsigned DIM>
double ImplicitCardiacMechanicsAssembler< DIM >::mNextTime [private]

Time to which the solver has been asked to solve to

Definition at line 81 of file ImplicitCardiacMechanicsAssembler.hpp.

Referenced by ImplicitCardiacMechanicsAssembler< DIM >::AssembleOnElement(), and ImplicitCardiacMechanicsAssembler< DIM >::Solve().

template<unsigned DIM>
double ImplicitCardiacMechanicsAssembler< DIM >::mOdeTimestep [private]

template<unsigned DIM>
bool ImplicitCardiacMechanicsAssembler< DIM >::mAllocatedMaterialLawMemory [private]

template<unsigned DIM>
unsigned ImplicitCardiacMechanicsAssembler< DIM >::mTotalQuadPoints [private]


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

Generated on Tue Aug 4 16:11:19 2009 for Chaste by  doxygen 1.5.5