#include <NonlinearElasticityAssembler.hpp>
Inherits AbstractNonlinearElasticityAssembler< DIM >.
Inherited by AbstractCardiacMechanicsAssembler< DIM >.
Public Member Functions | |
NonlinearElasticityAssembler (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) | |
NonlinearElasticityAssembler (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) | |
~NonlinearElasticityAssembler () | |
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 () |
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 | TestNonlinearElasticityAssembler |
Finite elasticity assembler. Solves static incompressible nonlinear elasticity problems with arbitrary material laws and a body force.
Uses quadratic-linear bases (for displacement and pressure), and is therefore outside the assembler 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 NonlinearElasticityAssembler.hpp.
NonlinearElasticityAssembler< DIM >::NonlinearElasticityAssembler | ( | 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 737 of file NonlinearElasticityAssembler.cpp.
References NonlinearElasticityAssembler< DIM >::Initialise().
NonlinearElasticityAssembler< DIM >::NonlinearElasticityAssembler | ( | 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 755 of file NonlinearElasticityAssembler.cpp.
References AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), and NonlinearElasticityAssembler< DIM >::Initialise().
NonlinearElasticityAssembler< DIM >::~NonlinearElasticityAssembler | ( | ) | [inline] |
Destructor frees memory for quadrature rules.
Definition at line 774 of file NonlinearElasticityAssembler.cpp.
References NonlinearElasticityAssembler< DIM >::mpBoundaryQuadratureRule, and NonlinearElasticityAssembler< DIM >::mpQuadratureRule.
void NonlinearElasticityAssembler< 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 646 of file NonlinearElasticityAssembler.cpp.
References PetscTools::GetNumProcs(), LinearSystem::GetOwnershipRange(), AbstractNonlinearElasticityAssembler< DIM >::mNumDofs, AbstractNonlinearElasticityAssembler< DIM >::mpLinearSystem, AbstractNonlinearElasticityAssembler< DIM >::mpPreconditionMatrixLinearSystem, NonlinearElasticityAssembler< DIM >::mpQuadMesh, and LinearSystem::rGetLhsMatrix().
Referenced by NonlinearElasticityAssembler< DIM >::Initialise().
void NonlinearElasticityAssembler< 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 514 of file NonlinearElasticityAssembler.cpp.
References AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), GaussianQuadratureRule< ELEMENT_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEMENT_DIM >::GetWeight(), NonlinearElasticityAssembler< DIM >::mpBoundaryQuadratureRule, NonlinearElasticityAssembler< DIM >::mpQuadMesh, AbstractNonlinearElasticityAssembler< DIM >::mpTractionBoundaryConditionFunction, AbstractNonlinearElasticityAssembler< DIM >::mUsingTractionBoundaryConditionFunction, NonlinearElasticityAssembler< DIM >::NUM_NODES_PER_BOUNDARY_ELEMENT, and GaussianQuadratureRule< ELEMENT_DIM >::rGetQuadPoint().
Referenced by NonlinearElasticityAssembler< DIM >::AssembleSystem().
void NonlinearElasticityAssembler< 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. |
Reimplemented in AbstractCardiacMechanicsAssembler< DIM >.
Definition at line 191 of file NonlinearElasticityAssembler.cpp.
References AbstractIncompressibleMaterialLaw< DIM >::ComputeStressAndStressDerivative(), AbstractNonlinearElasticityAssembler< DIM >::dTdE, AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), GaussianQuadratureRule< ELEMENT_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEMENT_DIM >::GetWeight(), AbstractNonlinearElasticityAssembler< DIM >::mBodyForce, AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, AbstractNonlinearElasticityAssembler< DIM >::mDensity, AbstractNonlinearElasticityAssembler< DIM >::mMaterialLaws, AbstractNonlinearElasticityAssembler< DIM >::mpBodyForceFunction, NonlinearElasticityAssembler< DIM >::mpQuadMesh, NonlinearElasticityAssembler< DIM >::mpQuadratureRule, AbstractNonlinearElasticityAssembler< DIM >::mUsingBodyForceFunction, NonlinearElasticityAssembler< DIM >::NUM_NODES_PER_ELEMENT, NonlinearElasticityAssembler< DIM >::NUM_VERTICES_PER_ELEMENT, GaussianQuadratureRule< ELEMENT_DIM >::rGetQuadPoint(), FourthOrderTensor< DIM >::SetAsProduct(), and NonlinearElasticityAssembler< DIM >::STENCIL_SIZE.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleSystem().
void NonlinearElasticityAssembler< 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 AbstractNonlinearElasticityAssembler< DIM >.
Definition at line 47 of file NonlinearElasticityAssembler.cpp.
References LinearSystem::AddLhsMultipleValues(), LinearSystem::AddRhsMultipleValues(), AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), LinearSystem::AssembleFinalLhsMatrix(), LinearSystem::AssembleIntermediateLhsMatrix(), NonlinearElasticityAssembler< DIM >::AssembleOnBoundaryElement(), NonlinearElasticityAssembler< DIM >::AssembleOnElement(), LinearSystem::AssembleRhsVector(), NonlinearElasticityAssembler< DIM >::BOUNDARY_STENCIL_SIZE, PetscTools::GetMyRank(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetOwnership(), NonlinearElasticityAssembler< DIM >::mBoundaryElements, AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, AbstractNonlinearElasticityAssembler< DIM >::mNumDofs, AbstractNonlinearElasticityAssembler< DIM >::mpLinearSystem, AbstractNonlinearElasticityAssembler< DIM >::mpPreconditionMatrixLinearSystem, NonlinearElasticityAssembler< DIM >::mpQuadMesh, AbstractNonlinearElasticityAssembler< DIM >::mSurfaceTractions, NonlinearElasticityAssembler< DIM >::NUM_NODES_PER_BOUNDARY_ELEMENT, NonlinearElasticityAssembler< DIM >::NUM_NODES_PER_ELEMENT, NonlinearElasticityAssembler< DIM >::NUM_VERTICES_PER_ELEMENT, NonlinearElasticityAssembler< DIM >::STENCIL_SIZE, LinearSystem::ZeroLhsMatrix(), and LinearSystem::ZeroRhsVector().
Referenced by ImplicitCardiacMechanicsAssembler< DIM >::Solve(), and ExplicitCardiacMechanicsAssembler< DIM >::Solve().
void NonlinearElasticityAssembler< 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 AbstractNonlinearElasticityAssembler< DIM >.
Definition at line 578 of file NonlinearElasticityAssembler.cpp.
References AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, AbstractNonlinearElasticityAssembler< DIM >::mMaterialLaws, AbstractNonlinearElasticityAssembler< DIM >::mNumDofs, NonlinearElasticityAssembler< DIM >::mpQuadMesh, and NonlinearElasticityAssembler< DIM >::NUM_VERTICES_PER_ELEMENT.
Referenced by NonlinearElasticityAssembler< DIM >::Initialise().
void NonlinearElasticityAssembler< DIM >::Initialise | ( | std::vector< c_vector< double, DIM > > * | pFixedNodeLocations | ) | [inline, protected] |
Initialise the assembler.
pFixedNodeLocations |
Definition at line 606 of file NonlinearElasticityAssembler.cpp.
References NonlinearElasticityAssembler< DIM >::AllocateMatrixMemory(), NonlinearElasticityAssembler< DIM >::FormInitialGuess(), AbstractNonlinearElasticityAssembler< DIM >::mFixedNodeDisplacements, AbstractNonlinearElasticityAssembler< DIM >::mFixedNodes, NonlinearElasticityAssembler< DIM >::mpBoundaryQuadratureRule, NonlinearElasticityAssembler< DIM >::mpQuadMesh, and NonlinearElasticityAssembler< DIM >::mpQuadratureRule.
Referenced by NonlinearElasticityAssembler< DIM >::NonlinearElasticityAssembler().
std::vector< c_vector< double, DIM > > & NonlinearElasticityAssembler< DIM >::rGetDeformedPosition | ( | ) | [inline, virtual] |
Get the deformed position. Note returnvalue[i](j) = x_j for node i.
Implements AbstractNonlinearElasticityAssembler< DIM >.
Definition at line 814 of file NonlinearElasticityAssembler.cpp.
References AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, AbstractNonlinearElasticityAssembler< DIM >::mDeformedPosition, and NonlinearElasticityAssembler< DIM >::mpQuadMesh.
std::vector< double > & NonlinearElasticityAssembler< DIM >::rGetPressures | ( | ) | [inline] |
Get pressures.
Definition at line 801 of file NonlinearElasticityAssembler.cpp.
References AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, NonlinearElasticityAssembler< DIM >::mpQuadMesh, and AbstractNonlinearElasticityAssembler< DIM >::mPressures.
void NonlinearElasticityAssembler< 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 791 of file NonlinearElasticityAssembler.cpp.
References NonlinearElasticityAssembler< DIM >::mBoundaryElements, AbstractNonlinearElasticityAssembler< DIM >::mpTractionBoundaryConditionFunction, and AbstractNonlinearElasticityAssembler< DIM >::mUsingTractionBoundaryConditionFunction.
void NonlinearElasticityAssembler< 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 781 of file NonlinearElasticityAssembler.cpp.
References NonlinearElasticityAssembler< DIM >::mBoundaryElements, and AbstractNonlinearElasticityAssembler< DIM >::mSurfaceTractions.
const size_t NonlinearElasticityAssembler< DIM >::BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + DIM [static, protected] |
Boundary stencil size
Definition at line 76 of file NonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleSystem().
std::vector<BoundaryElement<DIM-1,DIM>*> NonlinearElasticityAssembler< DIM >::mBoundaryElements [protected] |
Boundary elements with (non-zero) surface tractions defined on them
Definition at line 85 of file NonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleSystem(), NonlinearElasticityAssembler< DIM >::SetFunctionalTractionBoundaryCondition(), and NonlinearElasticityAssembler< DIM >::SetSurfaceTractionBoundaryConditions().
GaussianQuadratureRule<DIM-1>* NonlinearElasticityAssembler< DIM >::mpBoundaryQuadratureRule [protected] |
Boundary Gaussian quadrature rule
Definition at line 91 of file NonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnBoundaryElement(), NonlinearElasticityAssembler< DIM >::Initialise(), and NonlinearElasticityAssembler< DIM >::~NonlinearElasticityAssembler().
QuadraticMesh<DIM>* NonlinearElasticityAssembler< 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 NonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AllocateMatrixMemory(), NonlinearElasticityAssembler< DIM >::AssembleOnBoundaryElement(), NonlinearElasticityAssembler< DIM >::AssembleOnElement(), AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), NonlinearElasticityAssembler< DIM >::AssembleSystem(), NonlinearElasticityAssembler< DIM >::FormInitialGuess(), NonlinearElasticityAssembler< DIM >::Initialise(), NonlinearElasticityAssembler< DIM >::rGetDeformedPosition(), NonlinearElasticityAssembler< DIM >::rGetPressures(), and AbstractCardiacMechanicsAssembler< DIM >::SetVariableFibreSheetDirections().
GaussianQuadratureRule<DIM>* NonlinearElasticityAssembler< DIM >::mpQuadratureRule [protected] |
Gaussian quadrature rule
Definition at line 88 of file NonlinearElasticityAssembler.hpp.
Referenced by AbstractCardiacMechanicsAssembler< DIM >::AbstractCardiacMechanicsAssembler(), NonlinearElasticityAssembler< DIM >::AssembleOnElement(), AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), AbstractCardiacMechanicsAssembler< DIM >::GetQuadratureRule(), NonlinearElasticityAssembler< DIM >::Initialise(), and NonlinearElasticityAssembler< DIM >::~NonlinearElasticityAssembler().
const size_t NonlinearElasticityAssembler< DIM >::NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2 [static, protected] |
Number of nodes per boundary element
Definition at line 74 of file NonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnBoundaryElement(), and NonlinearElasticityAssembler< DIM >::AssembleSystem().
const size_t NonlinearElasticityAssembler< DIM >::NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2 [static, protected] |
Number of nodes per element
Reimplemented in AbstractCardiacMechanicsAssembler< DIM >.
Definition at line 70 of file NonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnElement(), and NonlinearElasticityAssembler< DIM >::AssembleSystem().
const size_t NonlinearElasticityAssembler< DIM >::NUM_VERTICES_PER_ELEMENT = DIM+1 [static, protected] |
Number of vertices per element
Reimplemented in AbstractCardiacMechanicsAssembler< DIM >.
Definition at line 68 of file NonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnElement(), NonlinearElasticityAssembler< DIM >::AssembleSystem(), and NonlinearElasticityAssembler< DIM >::FormInitialGuess().
const size_t NonlinearElasticityAssembler< DIM >::STENCIL_SIZE = DIM*NUM_NODES_PER_ELEMENT + NUM_VERTICES_PER_ELEMENT [static, protected] |
Stencil size
Reimplemented in AbstractCardiacMechanicsAssembler< DIM >.
Definition at line 72 of file NonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnElement(), and NonlinearElasticityAssembler< DIM >::AssembleSystem().