NonlinearElasticitySolver< DIM > Class Template Reference

#include <NonlinearElasticitySolver.hpp>

Inheritance diagram for NonlinearElasticitySolver< DIM >:

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

Collaboration graph
[legend]

List of all members.

Public Member Functions

 NonlinearElasticitySolver (QuadraticMesh< DIM > *pQuadMesh, AbstractMaterialLaw< 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< AbstractMaterialLaw< 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 ()
std::vector< double > & rGetPressures ()

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)

Protected Attributes

std::vector
< AbstractIncompressibleMaterialLaw
< DIM > * > 
mMaterialLaws
std::vector< double > mPressures

Static Protected Attributes

static const size_t NUM_NODES_PER_ELEMENT = AbstractNonlinearElasticitySolver<DIM>::NUM_NODES_PER_ELEMENT
static const size_t NUM_VERTICES_PER_ELEMENT = AbstractNonlinearElasticitySolver<DIM>::NUM_VERTICES_PER_ELEMENT
static const size_t NUM_NODES_PER_BOUNDARY_ELEMENT = AbstractNonlinearElasticitySolver<DIM>::NUM_NODES_PER_BOUNDARY_ELEMENT
static const size_t STENCIL_SIZE = DIM*NUM_NODES_PER_ELEMENT + NUM_VERTICES_PER_ELEMENT
static const size_t BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + DIM

Friends

class TestNonlinearElasticitySolver
class TestNonlinearElasticityAdjointSolver
class AdaptiveNonlinearElasticityProblem


Detailed Description

template<size_t DIM>
class NonlinearElasticitySolver< DIM >

Finite elasticity solver. Solves static *incompressible* nonlinear elasticity problems with arbitrary (incompressible) material laws and a body force.

Uses quadratic-linear bases (for displacement and pressure), and is therefore outside other assembler or solver hierarchy.

Definition at line 53 of file NonlinearElasticitySolver.hpp.


Constructor & Destructor Documentation

template<size_t DIM>
NonlinearElasticitySolver< DIM >::NonlinearElasticitySolver ( QuadraticMesh< DIM > *  pQuadMesh,
AbstractMaterialLaw< 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 for homogeneous problems.

Parameters:
pQuadMesh The quadratic mesh to solve on
pMaterialLaw A single material law to use on all elements
bodyForce The body force if constant. (If not constant, pass in a zero vector and call SetFunctionalBodyForce()
density The density (assumed constant)
outputDirectory The output directory
fixedNodes Which nodes are fixed in space (the displacement is assumed to be zero unless the next parameter is given
pFixedNodeLocations Optional new locations of the fixed nodes.

Definition at line 676 of file NonlinearElasticitySolver.cpp.

References EXCEPTION, NonlinearElasticitySolver< DIM >::FormInitialGuess(), AbstractNonlinearElasticitySolver< DIM >::Initialise(), and NonlinearElasticitySolver< DIM >::mMaterialLaws.

template<size_t DIM>
NonlinearElasticitySolver< DIM >::NonlinearElasticitySolver ( QuadraticMesh< DIM > *  pQuadMesh,
std::vector< AbstractMaterialLaw< 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 for heterogeneous problems

Parameters:
pQuadMesh The quadratic mesh to solve on
rMaterialLaws Vector of material laws for each element
bodyForce The body force if constant. (If not constant, pass in a zero vector and call SetFunctionalBodyForce()
density The density (assumed constant)
outputDirectory The output directory
fixedNodes Which nodes are fixed in space (the displacement is assumed to be zero unless the next parameter is given
pFixedNodeLocations Optional new locations of the fixed nodes.

Definition at line 704 of file NonlinearElasticitySolver.cpp.

References EXCEPTION, NonlinearElasticitySolver< DIM >::FormInitialGuess(), AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::GetNumElements(), AbstractNonlinearElasticitySolver< DIM >::Initialise(), and NonlinearElasticitySolver< DIM >::mMaterialLaws.

template<size_t DIM>
NonlinearElasticitySolver< DIM >::~NonlinearElasticitySolver (  )  [inline]

Destructor.

Definition at line 736 of file NonlinearElasticitySolver.cpp.


Member Function Documentation

template<size_t DIM>
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].

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.

Definition at line 192 of file NonlinearElasticitySolver.cpp.

References QuadraticBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), AbstractNonlinearElasticitySolver< 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 >::mCurrentTime, AbstractNonlinearElasticitySolver< DIM >::mDensity, NonlinearElasticitySolver< DIM >::mMaterialLaws, AbstractNonlinearElasticitySolver< DIM >::mpBodyForceFunction, AbstractNonlinearElasticitySolver< DIM >::mpQuadMesh, AbstractNonlinearElasticitySolver< 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().

template<size_t DIM>
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.

Parameters:
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 575 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(), AbstractNonlinearElasticitySolver< DIM >::mCurrentTime, AbstractNonlinearElasticitySolver< DIM >::mpBoundaryQuadratureRule, AbstractNonlinearElasticitySolver< 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().

template<size_t DIM>
void NonlinearElasticitySolver< DIM >::FormInitialGuess (  )  [inline, protected]

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.

Definition at line 639 of file NonlinearElasticitySolver.cpp.

References AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, NonlinearElasticitySolver< DIM >::mMaterialLaws, AbstractNonlinearElasticitySolver< DIM >::mNumDofs, AbstractNonlinearElasticitySolver< DIM >::mpQuadMesh, and NonlinearElasticitySolver< DIM >::NUM_VERTICES_PER_ELEMENT.

Referenced by NonlinearElasticitySolver< DIM >::NonlinearElasticitySolver().

template<size_t DIM>
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).

Parameters:
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 NonlinearElasticitySolver< DIM >::AssembleOnBoundaryElement(), NonlinearElasticitySolver< DIM >::AssembleOnElement(), NonlinearElasticitySolver< DIM >::BOUNDARY_STENCIL_SIZE, AbstractNonlinearElasticitySolver< DIM >::FinishAssembleSystem(), PetscTools::GetMyRank(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetOwnership(), AbstractNonlinearElasticitySolver< DIM >::mBoundaryElements, AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< DIM >::mJacobianMatrix, AbstractNonlinearElasticitySolver< DIM >::mNumDofs, AbstractNonlinearElasticitySolver< DIM >::mpQuadMesh, AbstractNonlinearElasticitySolver< DIM >::mPreconditionMatrix, AbstractNonlinearElasticitySolver< DIM >::mResidualVector, 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, PetscMatTools::Zero(), and PetscVecTools::Zero().

template<size_t DIM>
std::vector< double > & NonlinearElasticitySolver< DIM >::rGetPressures (  )  [inline]


Member Data Documentation

template<size_t DIM>
const size_t NonlinearElasticitySolver< DIM >::NUM_NODES_PER_ELEMENT = AbstractNonlinearElasticitySolver<DIM>::NUM_NODES_PER_ELEMENT [static, protected]

template<size_t DIM>
const size_t NonlinearElasticitySolver< DIM >::NUM_VERTICES_PER_ELEMENT = AbstractNonlinearElasticitySolver<DIM>::NUM_VERTICES_PER_ELEMENT [static, protected]

template<size_t DIM>
const size_t NonlinearElasticitySolver< DIM >::STENCIL_SIZE = DIM*NUM_NODES_PER_ELEMENT + NUM_VERTICES_PER_ELEMENT [static, protected]

Stencil size - number of unknowns per element (DIM*NUM_NODES_PER_ELEMENT displacement unknowns, NUM_VERTICES_PER_ELEMENT pressure unknowns

Definition at line 69 of file NonlinearElasticitySolver.hpp.

Referenced by NonlinearElasticitySolver< DIM >::AssembleOnElement(), and NonlinearElasticitySolver< DIM >::AssembleSystem().

template<size_t DIM>
const size_t NonlinearElasticitySolver< DIM >::BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + DIM [static, protected]

Boundary stencil size

Definition at line 71 of file NonlinearElasticitySolver.hpp.

Referenced by NonlinearElasticitySolver< DIM >::AssembleSystem().

template<size_t DIM>
std::vector<AbstractIncompressibleMaterialLaw<DIM>*> NonlinearElasticitySolver< DIM >::mMaterialLaws [protected]

The material laws for each element. This will either be of size 1 (same material law for all elements, ie homogeneous), or size num_elem.

Definition at line 79 of file NonlinearElasticitySolver.hpp.

Referenced by NonlinearElasticitySolver< DIM >::AssembleOnElement(), NonlinearElasticitySolver< DIM >::FormInitialGuess(), and NonlinearElasticitySolver< DIM >::NonlinearElasticitySolver().

template<size_t DIM>
std::vector<double> NonlinearElasticitySolver< DIM >::mPressures [protected]

The solution pressures. mPressures[i] = pressure at node i (ie vertex i).

Definition at line 84 of file NonlinearElasticitySolver.hpp.

Referenced by NonlinearElasticitySolver< DIM >::rGetPressures().


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

Generated on Tue May 31 14:33:58 2011 for Chaste by  doxygen 1.5.5