#include <NonlinearElasticitySolver.hpp>
Public Member Functions | |
NonlinearElasticitySolver (QuadraticMesh< DIM > *pQuadMesh, AbstractIncompressibleMaterialLaw< DIM > *pMaterialLaw, c_vector< double, DIM > bodyForce, double density, std::string outputDirectory, std::vector< unsigned > &fixedNodes, std::vector< c_vector< double, DIM > > *pFixedNodeLocations=NULL) | |
NonlinearElasticitySolver (QuadraticMesh< DIM > *pQuadMesh, std::vector< AbstractIncompressibleMaterialLaw< DIM > * > &rMaterialLaws, c_vector< double, DIM > bodyForce, double density, std::string outputDirectory, std::vector< unsigned > &fixedNodes, std::vector< c_vector< double, DIM > > *pFixedNodeLocations=NULL) | |
~NonlinearElasticitySolver () | |
void | SetSurfaceTractionBoundaryConditions (std::vector< BoundaryElement< DIM-1, DIM > * > &rBoundaryElements, std::vector< c_vector< double, DIM > > &rSurfaceTractions) |
void | SetFunctionalTractionBoundaryCondition (std::vector< BoundaryElement< DIM-1, DIM > * > rBoundaryElements, c_vector< double, DIM >(*pFunction)(c_vector< double, DIM > &)) |
std::vector< double > & | rGetPressures () |
std::vector< c_vector< double, DIM > > & | rGetDeformedPosition () |
Protected Member Functions | |
virtual 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 | AssembleOnBoundaryElement (BoundaryElement< DIM-1, DIM > &rBoundaryElement, c_matrix< double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE > &rAelem, c_vector< double, BOUNDARY_STENCIL_SIZE > &rBelem, c_vector< double, DIM > &rTraction, bool assembleResidual, bool assembleJacobian) |
void | FormInitialGuess () |
void | AssembleSystem (bool assembleResidual, bool assembleJacobian) |
void | Initialise (std::vector< c_vector< double, DIM > > *pFixedNodeLocations) |
void | AllocateMatrixMemory () |
virtual 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) |
Protected Attributes | |
QuadraticMesh< DIM > * | mpQuadMesh |
std::vector< BoundaryElement < DIM-1, DIM > * > | mBoundaryElements |
GaussianQuadratureRule< DIM > * | mpQuadratureRule |
GaussianQuadratureRule< DIM-1 > * | mpBoundaryQuadratureRule |
Static Protected Attributes | |
static const size_t | NUM_VERTICES_PER_ELEMENT = DIM+1 |
static const size_t | NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2 |
static const size_t | STENCIL_SIZE = DIM*NUM_NODES_PER_ELEMENT + NUM_VERTICES_PER_ELEMENT |
static const size_t | NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2 |
static const size_t | BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + DIM |
Friends | |
class | TestNonlinearElasticitySolver |
Uses quadratic-linear bases (for displacement and pressure), and is therefore outside other assember or solver hierachy.
Currently only works with fixed nodes BCs (ie zerodisplacement) and zero-surface tractions on the rest of the boundary.
Definition at line 61 of file NonlinearElasticitySolver.hpp.
NonlinearElasticitySolver< DIM >::NonlinearElasticitySolver | ( | QuadraticMesh< DIM > * | pQuadMesh, | |
AbstractIncompressibleMaterialLaw< DIM > * | pMaterialLaw, | |||
c_vector< double, DIM > | bodyForce, | |||
double | density, | |||
std::string | outputDirectory, | |||
std::vector< unsigned > & | fixedNodes, | |||
std::vector< c_vector< double, DIM > > * | pFixedNodeLocations = NULL | |||
) | [inline] |
Constructor taking in mesh, material law (assuming homogeniety at the moment) body force, density, the fixed nodes (all the fixed nodes, including non-vertices), and the output directory.
pQuadMesh | ||
pMaterialLaw | ||
bodyForce | ||
density | ||
outputDirectory | ||
fixedNodes | ||
pFixedNodeLocations | (defaults to NULL) |
Definition at line 810 of file NonlinearElasticitySolver.cpp.
References NonlinearElasticitySolver< DIM >::Initialise().
NonlinearElasticitySolver< DIM >::NonlinearElasticitySolver | ( | QuadraticMesh< DIM > * | pQuadMesh, | |
std::vector< AbstractIncompressibleMaterialLaw< DIM > * > & | rMaterialLaws, | |||
c_vector< double, DIM > | bodyForce, | |||
double | density, | |||
std::string | outputDirectory, | |||
std::vector< unsigned > & | fixedNodes, | |||
std::vector< c_vector< double, DIM > > * | pFixedNodeLocations = NULL | |||
) | [inline] |
Variant constructor taking a vector of material laws.
pQuadMesh | ||
rMaterialLaws | ||
bodyForce | ||
density | ||
outputDirectory | ||
fixedNodes | ||
pFixedNodeLocations | (defaults to NULL) |
Definition at line 828 of file NonlinearElasticitySolver.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), and NonlinearElasticitySolver< DIM >::Initialise().
NonlinearElasticitySolver< DIM >::~NonlinearElasticitySolver | ( | ) | [inline] |
Destructor frees memory for quadrature rules.
Definition at line 847 of file NonlinearElasticitySolver.cpp.
References NonlinearElasticitySolver< DIM >::mpBoundaryQuadratureRule, and NonlinearElasticitySolver< DIM >::mpQuadratureRule.
void NonlinearElasticitySolver< 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] |
Assemble residual or jacobian on an element, using the current solution stored in mCurrrentSolution. The ordering assumed is (in 2d) rBElem = [u0 v0 u1 v1 .. u5 v5 p0 p1 p2].
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. |
Definition at line 204 of file NonlinearElasticitySolver.cpp.
References QuadraticBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), NonlinearElasticitySolver< DIM >::ComputeStressAndStressDerivative(), QuadraticBasisFunction< ELEMENT_DIM >::ComputeTransformedBasisFunctionDerivatives(), Determinant(), AbstractNonlinearElasticitySolver< DIM >::dTdE, AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), GaussianQuadratureRule< ELEMENT_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEMENT_DIM >::GetWeight(), Inverse(), AbstractNonlinearElasticitySolver< DIM >::mBodyForce, AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< DIM >::mDensity, AbstractNonlinearElasticitySolver< DIM >::mMaterialLaws, AbstractNonlinearElasticitySolver< DIM >::mpBodyForceFunction, NonlinearElasticitySolver< DIM >::mpQuadMesh, NonlinearElasticitySolver< DIM >::mpQuadratureRule, AbstractNonlinearElasticitySolver< DIM >::mUsingBodyForceFunction, NonlinearElasticitySolver< DIM >::NUM_NODES_PER_ELEMENT, NonlinearElasticitySolver< DIM >::NUM_VERTICES_PER_ELEMENT, GaussianQuadratureRule< ELEMENT_DIM >::rGetQuadPoint(), and NonlinearElasticitySolver< DIM >::STENCIL_SIZE.
Referenced by NonlinearElasticitySolver< DIM >::AssembleSystem().
void NonlinearElasticitySolver< DIM >::AssembleOnBoundaryElement | ( | BoundaryElement< DIM-1, DIM > & | rBoundaryElement, | |
c_matrix< double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE > & | rAelem, | |||
c_vector< double, BOUNDARY_STENCIL_SIZE > & | rBelem, | |||
c_vector< double, DIM > & | rTraction, | |||
bool | assembleResidual, | |||
bool | assembleJacobian | |||
) | [inline, protected, virtual] |
Compute the term from the surface integral of s*phi, where s is a specified non-zero surface traction (ie Neumann boundary condition) to be added to the Rhs vector.
rBoundaryElement | ||
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. | |
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. | |
rTraction | ||
assembleResidual | A bool stating whether to assemble the residual vector. | |
assembleJacobian | A bool stating whether to assemble the Jacobian matrix. |
Definition at line 590 of file NonlinearElasticitySolver.cpp.
References QuadraticBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), GaussianQuadratureRule< ELEMENT_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEMENT_DIM >::GetWeight(), NonlinearElasticitySolver< DIM >::mpBoundaryQuadratureRule, NonlinearElasticitySolver< DIM >::mpQuadMesh, AbstractNonlinearElasticitySolver< DIM >::mpTractionBoundaryConditionFunction, AbstractNonlinearElasticitySolver< DIM >::mUsingTractionBoundaryConditionFunction, NonlinearElasticitySolver< DIM >::NUM_NODES_PER_BOUNDARY_ELEMENT, and GaussianQuadratureRule< ELEMENT_DIM >::rGetQuadPoint().
Referenced by NonlinearElasticitySolver< DIM >::AssembleSystem().
void NonlinearElasticitySolver< DIM >::FormInitialGuess | ( | ) | [inline, protected, virtual] |
Set up the current guess to be the solution given no displacement. The current solution (in 2d) is order as [u1 v1 u2 v2 ... uN vN p1 p2 .. pM] (where there are N total nodes and M vertices) so the initial guess is [0 0 0 0 ... 0 0 p1 p2 .. pM] where p_i are such that T is zero (depends on material law).
In a homogeneous problem, all p_i are the same. In a heterogeneous problem, p for a given vertex is the zero-strain-pressure for ONE of the elements containing that vertex (which element containing the vertex is reached LAST). In this case the initial guess will be close but not exactly the solution given zero body force.
Implements AbstractNonlinearElasticitySolver< DIM >.
Definition at line 654 of file NonlinearElasticitySolver.cpp.
References AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< DIM >::mMaterialLaws, AbstractNonlinearElasticitySolver< DIM >::mNumDofs, NonlinearElasticitySolver< DIM >::mpQuadMesh, and NonlinearElasticitySolver< DIM >::NUM_VERTICES_PER_ELEMENT.
Referenced by NonlinearElasticitySolver< DIM >::Initialise().
void NonlinearElasticitySolver< DIM >::AssembleSystem | ( | bool | assembleResidual, | |
bool | assembleJacobian | |||
) | [inline, protected, virtual] |
Assemble the residual vector (using the current solution stored in mCurrentSolution, output going to mpLinearSystem->rGetRhsVector), or Jacobian matrix (using the current solution stored in mCurrentSolution, output going to mpLinearSystem->rGetLhsMatrix).
assembleResidual | A bool stating whether to assemble the residual vector. | |
assembleJacobian | A bool stating whether to assemble the Jacobian matrix. |
Implements AbstractNonlinearElasticitySolver< DIM >.
Definition at line 48 of file NonlinearElasticitySolver.cpp.
References LinearSystem::AddLhsMultipleValues(), LinearSystem::AddRhsMultipleValues(), AbstractNonlinearElasticitySolver< DIM >::ApplyBoundaryConditions(), LinearSystem::AssembleFinalLhsMatrix(), LinearSystem::AssembleIntermediateLhsMatrix(), NonlinearElasticitySolver< DIM >::AssembleOnBoundaryElement(), NonlinearElasticitySolver< DIM >::AssembleOnElement(), LinearSystem::AssembleRhsVector(), NonlinearElasticitySolver< DIM >::BOUNDARY_STENCIL_SIZE, PetscTools::GetMyRank(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetOwnership(), NonlinearElasticitySolver< DIM >::mBoundaryElements, AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< DIM >::mNumDofs, AbstractNonlinearElasticitySolver< DIM >::mpLinearSystem, AbstractNonlinearElasticitySolver< DIM >::mpPreconditionMatrixLinearSystem, NonlinearElasticitySolver< DIM >::mpQuadMesh, AbstractNonlinearElasticitySolver< DIM >::mSurfaceTractions, NonlinearElasticitySolver< DIM >::NUM_NODES_PER_BOUNDARY_ELEMENT, NonlinearElasticitySolver< DIM >::NUM_NODES_PER_ELEMENT, NonlinearElasticitySolver< DIM >::NUM_VERTICES_PER_ELEMENT, NonlinearElasticitySolver< DIM >::STENCIL_SIZE, LinearSystem::ZeroLhsMatrix(), and LinearSystem::ZeroRhsVector().
Referenced by ImplicitCardiacMechanicsSolver< DIM >::Solve(), and ExplicitCardiacMechanicsSolver< DIM >::Solve().
void NonlinearElasticitySolver< DIM >::Initialise | ( | std::vector< c_vector< double, DIM > > * | pFixedNodeLocations | ) | [inline, protected] |
Initialise the solver.
pFixedNodeLocations |
Definition at line 682 of file NonlinearElasticitySolver.cpp.
References NonlinearElasticitySolver< DIM >::AllocateMatrixMemory(), NonlinearElasticitySolver< DIM >::FormInitialGuess(), AbstractNonlinearElasticitySolver< DIM >::mFixedNodeDisplacements, AbstractNonlinearElasticitySolver< DIM >::mFixedNodes, NonlinearElasticitySolver< DIM >::mpBoundaryQuadratureRule, NonlinearElasticitySolver< DIM >::mpQuadMesh, and NonlinearElasticitySolver< DIM >::mpQuadratureRule.
Referenced by NonlinearElasticitySolver< DIM >::NonlinearElasticitySolver().
void NonlinearElasticitySolver< DIM >::AllocateMatrixMemory | ( | ) | [inline, protected] |
Allocates memory for the Jacobian and preconditioner matrices (larger number of non-zeros per row than with say linear problems)
Definition at line 722 of file NonlinearElasticitySolver.cpp.
References PetscTools::GetNumProcs(), LinearSystem::GetOwnershipRange(), AbstractNonlinearElasticitySolver< DIM >::mNumDofs, AbstractNonlinearElasticitySolver< DIM >::mpLinearSystem, AbstractNonlinearElasticitySolver< DIM >::mpPreconditionMatrixLinearSystem, NonlinearElasticitySolver< DIM >::mpQuadMesh, and LinearSystem::rGetLhsMatrix().
Referenced by NonlinearElasticitySolver< DIM >::Initialise().
void NonlinearElasticitySolver< 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] |
Simple (one-line function which just calls ComputeStressAndStressDerivative() on the material law, using C, inv(C), and p as the input and with rT and rDTdE as the output. Overloaded by other assemblers (eg cardiac mechanics) which need to add extra terms to the stress.
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 in AbstractCardiacMechanicsSolver< DIM >.
Definition at line 574 of file NonlinearElasticitySolver.cpp.
References AbstractIncompressibleMaterialLaw< DIM >::ComputeStressAndStressDerivative().
Referenced by NonlinearElasticitySolver< DIM >::AssembleOnElement().
void NonlinearElasticitySolver< DIM >::SetSurfaceTractionBoundaryConditions | ( | std::vector< BoundaryElement< DIM-1, DIM > * > & | rBoundaryElements, | |
std::vector< c_vector< double, DIM > > & | rSurfaceTractions | |||
) | [inline] |
Specify traction boundary conditions (if this is not called zero surface tractions are assumed. This method takes in a list of boundary elements and a corresponding list of surface tractions.
rBoundaryElements | ||
rSurfaceTractions |
Definition at line 885 of file NonlinearElasticitySolver.cpp.
References NonlinearElasticitySolver< DIM >::mBoundaryElements, and AbstractNonlinearElasticitySolver< DIM >::mSurfaceTractions.
void NonlinearElasticitySolver< DIM >::SetFunctionalTractionBoundaryCondition | ( | std::vector< BoundaryElement< DIM-1, DIM > * > | rBoundaryElements, | |
c_vector< double, DIM >(*)(c_vector< double, DIM > &) | pFunction | |||
) | [inline] |
Set a function which gives the surface traction as a function of X (undeformed position), together with the surface elements which make up the Neumann part of the boundary.
rBoundaryElements | ||
pFunction |
Definition at line 895 of file NonlinearElasticitySolver.cpp.
References NonlinearElasticitySolver< DIM >::mBoundaryElements, AbstractNonlinearElasticitySolver< DIM >::mpTractionBoundaryConditionFunction, and AbstractNonlinearElasticitySolver< DIM >::mUsingTractionBoundaryConditionFunction.
std::vector< double > & NonlinearElasticitySolver< DIM >::rGetPressures | ( | ) | [inline] |
Get pressures.
Definition at line 905 of file NonlinearElasticitySolver.cpp.
References AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, NonlinearElasticitySolver< DIM >::mpQuadMesh, and AbstractNonlinearElasticitySolver< DIM >::mPressures.
std::vector< c_vector< double, DIM > > & NonlinearElasticitySolver< DIM >::rGetDeformedPosition | ( | ) | [inline, virtual] |
Get the deformed position. Note returnvalue[i](j) = x_j for node i.
Implements AbstractNonlinearElasticitySolver< DIM >.
Definition at line 918 of file NonlinearElasticitySolver.cpp.
References AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< DIM >::mDeformedPosition, and NonlinearElasticitySolver< DIM >::mpQuadMesh.
const size_t NonlinearElasticitySolver< DIM >::NUM_VERTICES_PER_ELEMENT = DIM+1 [static, protected] |
Number of vertices per element
Reimplemented in AbstractCardiacMechanicsSolver< DIM >.
Definition at line 68 of file NonlinearElasticitySolver.hpp.
Referenced by NonlinearElasticitySolver< DIM >::AssembleOnElement(), NonlinearElasticitySolver< DIM >::AssembleSystem(), and NonlinearElasticitySolver< DIM >::FormInitialGuess().
const size_t NonlinearElasticitySolver< DIM >::NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2 [static, protected] |
Number of nodes per element
Reimplemented in AbstractCardiacMechanicsSolver< DIM >.
Definition at line 70 of file NonlinearElasticitySolver.hpp.
Referenced by NonlinearElasticitySolver< DIM >::AssembleOnElement(), and NonlinearElasticitySolver< DIM >::AssembleSystem().
const size_t NonlinearElasticitySolver< DIM >::STENCIL_SIZE = DIM*NUM_NODES_PER_ELEMENT + NUM_VERTICES_PER_ELEMENT [static, protected] |
Stencil size
Reimplemented in AbstractCardiacMechanicsSolver< DIM >.
Definition at line 72 of file NonlinearElasticitySolver.hpp.
Referenced by NonlinearElasticitySolver< DIM >::AssembleOnElement(), and NonlinearElasticitySolver< DIM >::AssembleSystem().
const size_t NonlinearElasticitySolver< DIM >::NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2 [static, protected] |
Number of nodes per boundary element
Definition at line 74 of file NonlinearElasticitySolver.hpp.
Referenced by NonlinearElasticitySolver< DIM >::AssembleOnBoundaryElement(), and NonlinearElasticitySolver< DIM >::AssembleSystem().
const size_t NonlinearElasticitySolver< DIM >::BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + DIM [static, protected] |
Boundary stencil size
Definition at line 76 of file NonlinearElasticitySolver.hpp.
Referenced by NonlinearElasticitySolver< DIM >::AssembleSystem().
QuadraticMesh<DIM>* NonlinearElasticitySolver< DIM >::mpQuadMesh [protected] |
The mesh to be solved on. Requires 6 nodes per triangle (or 10 per tetrahedron) as quadratic bases are used.
Definition at line 82 of file NonlinearElasticitySolver.hpp.
Referenced by NonlinearElasticitySolver< DIM >::AllocateMatrixMemory(), NonlinearElasticitySolver< DIM >::AssembleOnBoundaryElement(), NonlinearElasticitySolver< DIM >::AssembleOnElement(), NonlinearElasticitySolver< DIM >::AssembleSystem(), AbstractCardiacMechanicsSolver< DIM >::ComputeDeformationGradientAndStretchInEachElement(), NonlinearElasticitySolver< DIM >::FormInitialGuess(), NonlinearElasticitySolver< DIM >::Initialise(), NonlinearElasticitySolver< DIM >::rGetDeformedPosition(), NonlinearElasticitySolver< DIM >::rGetPressures(), and AbstractCardiacMechanicsSolver< DIM >::SetVariableFibreSheetDirections().
std::vector<BoundaryElement<DIM-1,DIM>*> NonlinearElasticitySolver< DIM >::mBoundaryElements [protected] |
Boundary elements with (non-zero) surface tractions defined on them
Definition at line 85 of file NonlinearElasticitySolver.hpp.
Referenced by NonlinearElasticitySolver< DIM >::AssembleSystem(), NonlinearElasticitySolver< DIM >::SetFunctionalTractionBoundaryCondition(), and NonlinearElasticitySolver< DIM >::SetSurfaceTractionBoundaryConditions().
GaussianQuadratureRule<DIM>* NonlinearElasticitySolver< DIM >::mpQuadratureRule [protected] |
Gaussian quadrature rule
Definition at line 88 of file NonlinearElasticitySolver.hpp.
Referenced by AbstractCardiacMechanicsSolver< DIM >::AbstractCardiacMechanicsSolver(), NonlinearElasticitySolver< DIM >::AssembleOnElement(), AbstractCardiacMechanicsSolver< DIM >::GetQuadratureRule(), NonlinearElasticitySolver< DIM >::Initialise(), and NonlinearElasticitySolver< DIM >::~NonlinearElasticitySolver().
GaussianQuadratureRule<DIM-1>* NonlinearElasticitySolver< DIM >::mpBoundaryQuadratureRule [protected] |
Boundary Gaussian quadrature rule
Definition at line 91 of file NonlinearElasticitySolver.hpp.
Referenced by NonlinearElasticitySolver< DIM >::AssembleOnBoundaryElement(), NonlinearElasticitySolver< DIM >::Initialise(), and NonlinearElasticitySolver< DIM >::~NonlinearElasticitySolver().