NonlinearElasticityAssembler< DIM > Class Template Reference

#include <NonlinearElasticityAssembler.hpp>

Inheritance diagram for NonlinearElasticityAssembler< DIM >:

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

Collaboration graph
[legend]

List of all members.

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)

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


Detailed Description

template<size_t DIM>
class NonlinearElasticityAssembler< DIM >

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 62 of file NonlinearElasticityAssembler.hpp.


Constructor & Destructor Documentation

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

Definition at line 638 of file NonlinearElasticityAssembler.cpp.

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

Definition at line 661 of file NonlinearElasticityAssembler.cpp.

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

Destructor frees memory for quadrature rules

Definition at line 685 of file NonlinearElasticityAssembler.cpp.


Member Function Documentation

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

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

Definition at line 508 of file NonlinearElasticityAssembler.cpp.

References GaussianQuadratureRule< ELEM_DIM >::GetNumQuadPoints(), GaussianQuadratureRule< ELEM_DIM >::GetWeight(), NonlinearElasticityAssembler< DIM >::mpQuadMesh, and GaussianQuadratureRule< ELEM_DIM >::rGetQuadPoint().

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

template<size_t DIM>
void NonlinearElasticityAssembler< 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 572 of file NonlinearElasticityAssembler.cpp.

References NonlinearElasticityAssembler< DIM >::mpQuadMesh.

template<size_t DIM>
void NonlinearElasticityAssembler< DIM >::AssembleSystem ( bool  assembleResidual,
bool  assembleJacobian 
) [inline, protected]

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).

Definition at line 48 of file NonlinearElasticityAssembler.cpp.

References NonlinearElasticityAssembler< DIM >::AssembleOnBoundaryElement(), NonlinearElasticityAssembler< DIM >::AssembleOnElement(), NonlinearElasticityAssembler< DIM >::mBoundaryElements, and NonlinearElasticityAssembler< DIM >::mpQuadMesh.

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

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

Definition at line 692 of file NonlinearElasticityAssembler.cpp.

References NonlinearElasticityAssembler< DIM >::mBoundaryElements.

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

Definition at line 702 of file NonlinearElasticityAssembler.cpp.

References NonlinearElasticityAssembler< DIM >::mBoundaryElements.

template<size_t DIM>
std::vector< c_vector< double, DIM > > & NonlinearElasticityAssembler< DIM >::rGetDeformedPosition (  )  [inline]

Get the deformed position. Note returnvalue[i](j) = x_j for node i

Definition at line 725 of file NonlinearElasticityAssembler.cpp.

References NonlinearElasticityAssembler< DIM >::mpQuadMesh.


Member Data Documentation

template<size_t DIM>
QuadraticMesh<DIM>* NonlinearElasticityAssembler< DIM >::mpQuadMesh [protected]

template<size_t DIM>
std::vector<BoundaryElement<DIM-1,DIM>*> NonlinearElasticityAssembler< DIM >::mBoundaryElements [protected]


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

Generated on Wed Mar 18 12:52:46 2009 for Chaste by  doxygen 1.5.5