AbstractNonlinearElasticityAssembler< DIM > Class Template Reference

#include <AbstractNonlinearElasticityAssembler.hpp>

Inherited by NonlinearElasticityAssembler< DIM >.

Collaboration diagram for AbstractNonlinearElasticityAssembler< DIM >:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 AbstractNonlinearElasticityAssembler (unsigned numDofs, AbstractIncompressibleMaterialLaw< DIM > *pMaterialLaw, c_vector< double, DIM > bodyForce, double density, std::string outputDirectory, std::vector< unsigned > &fixedNodes)
 AbstractNonlinearElasticityAssembler (unsigned numDofs, std::vector< AbstractIncompressibleMaterialLaw< DIM > * > &rMaterialLaws, c_vector< double, DIM > bodyForce, double density, std::string outputDirectory, std::vector< unsigned > &fixedNodes)
virtual ~AbstractNonlinearElasticityAssembler ()
void Solve (double tol=-1.0, unsigned offset=0, unsigned maxNumNewtonIterations=INT_MAX, bool quitIfNoConvergence=true)
void WriteOutput (unsigned counter)
unsigned GetNumNewtonIterations ()
void SetFunctionalBodyForce (c_vector< double, DIM >(*pFunction)(c_vector< double, DIM > &))
void SetWriteOutput (bool writeOutput=true)

Protected Member Functions

virtual void FormInitialGuess ()=0
virtual void AssembleSystem (bool assembleResidual, bool assembleJacobian)=0
virtual std::vector< c_vector
< double, DIM > > & 
rGetDeformedPosition ()=0
void ApplyBoundaryConditions (bool applyToMatrix)
double ComputeResidualAndGetNorm ()
double CalculateResidualNorm ()
void VectorSum (std::vector< double > &rX, ReplicatableVector &rY, double a, std::vector< double > &rZ)
void PrintLineSearchResult (double s, double residNorm)
double TakeNewtonStep ()
double UpdateSolutionUsingLineSearch (Vec solution)
virtual void PostNewtonStep (unsigned counter, double normResidual)

Protected Attributes

unsigned mNumDofs
std::vector
< AbstractIncompressibleMaterialLaw
< DIM > * > 
mMaterialLaws
LinearSystemmpLinearSystem
LinearSystemmpPreconditionMatrixLinearSystem
c_vector< double, DIM > mBodyForce
double mDensity
std::string mOutputDirectory
std::vector< unsignedmFixedNodes
std::vector< c_vector< double,
DIM > > 
mFixedNodeDisplacements
bool mWriteOutput
std::vector< doublemCurrentSolution
FourthOrderTensor< DIM > dTdE
unsigned mNumNewtonIterations
std::vector< c_vector< double,
DIM > > 
mDeformedPosition
std::vector< doublemPressures
std::vector< c_vector< double,
DIM > > 
mSurfaceTractions
c_vector< double, DIM >(* mpBodyForceFunction )(c_vector< double, DIM > &)
c_vector< double, DIM >(* mpTractionBoundaryConditionFunction )(c_vector< double, DIM > &)
bool mUsingBodyForceFunction
bool mUsingTractionBoundaryConditionFunction

Static Protected Attributes

static const double MAX_NEWTON_ABS_TOL = 1e-7
static const double MIN_NEWTON_ABS_TOL = 1e-10
static const double NEWTON_REL_TOL = 1e-4

Detailed Description

template<unsigned DIM>
class AbstractNonlinearElasticityAssembler< DIM >

Abstract nonlinear elasticity assembler.

Definition at line 54 of file AbstractNonlinearElasticityAssembler.hpp.


Constructor & Destructor Documentation

template<unsigned DIM>
AbstractNonlinearElasticityAssembler< DIM >::AbstractNonlinearElasticityAssembler ( unsigned  numDofs,
AbstractIncompressibleMaterialLaw< DIM > *  pMaterialLaw,
c_vector< double, DIM >  bodyForce,
double  density,
std::string  outputDirectory,
std::vector< unsigned > &  fixedNodes 
) [inline]

Constructor.

Parameters:
numDofs 
pMaterialLaw 
bodyForce 
density 
outputDirectory 
fixedNodes 

Definition at line 699 of file AbstractNonlinearElasticityAssembler.hpp.

References AbstractNonlinearElasticityAssembler< DIM >::mMaterialLaws, AbstractNonlinearElasticityAssembler< DIM >::mOutputDirectory, and AbstractNonlinearElasticityAssembler< DIM >::mWriteOutput.

template<unsigned DIM>
AbstractNonlinearElasticityAssembler< DIM >::AbstractNonlinearElasticityAssembler ( unsigned  numDofs,
std::vector< AbstractIncompressibleMaterialLaw< DIM > * > &  rMaterialLaws,
c_vector< double, DIM >  bodyForce,
double  density,
std::string  outputDirectory,
std::vector< unsigned > &  fixedNodes 
) [inline]

Variant constructor taking a vector of material laws.

Parameters:
numDofs 
rMaterialLaws 
bodyForce 
density 
outputDirectory 
fixedNodes 

Definition at line 725 of file AbstractNonlinearElasticityAssembler.hpp.

References AbstractNonlinearElasticityAssembler< DIM >::mMaterialLaws, AbstractNonlinearElasticityAssembler< DIM >::mOutputDirectory, and AbstractNonlinearElasticityAssembler< DIM >::mWriteOutput.

template<unsigned DIM>
AbstractNonlinearElasticityAssembler< DIM >::~AbstractNonlinearElasticityAssembler (  )  [inline, virtual]

Member Function Documentation

template<unsigned DIM>
void AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions ( bool  applyToMatrix  )  [inline, protected]
template<unsigned DIM>
virtual void AbstractNonlinearElasticityAssembler< DIM >::AssembleSystem ( bool  assembleResidual,
bool  assembleJacobian 
) [protected, pure 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). Must be overridden in concrete derived classes.

Parameters:
assembleResidual A bool stating whether to assemble the residual vector.
assembleJacobian A bool stating whether to assemble the Jacobian matrix.

Implemented in NonlinearElasticityAssembler< DIM >.

Referenced by AbstractNonlinearElasticityAssembler< DIM >::ComputeResidualAndGetNorm(), and AbstractNonlinearElasticityAssembler< DIM >::TakeNewtonStep().

template<unsigned DIM>
double AbstractNonlinearElasticityAssembler< DIM >::CalculateResidualNorm (  )  [inline, protected]
template<unsigned DIM>
double AbstractNonlinearElasticityAssembler< DIM >::ComputeResidualAndGetNorm (  )  [inline, protected]

Set up the residual vector (using the current solution), and get its scaled norm (Calculate |r|_2 / length(r), where r is residual vector)

Definition at line 388 of file AbstractNonlinearElasticityAssembler.hpp.

References AbstractNonlinearElasticityAssembler< DIM >::AssembleSystem(), and AbstractNonlinearElasticityAssembler< DIM >::CalculateResidualNorm().

Referenced by AbstractNonlinearElasticityAssembler< DIM >::Solve(), and AbstractNonlinearElasticityAssembler< DIM >::UpdateSolutionUsingLineSearch().

template<unsigned DIM>
virtual void AbstractNonlinearElasticityAssembler< DIM >::FormInitialGuess (  )  [protected, pure virtual]

Set up the current guess to be the solution given no displacement. Must be overridden in concrete derived classes.

Implemented in NonlinearElasticityAssembler< DIM >.

template<unsigned DIM>
unsigned AbstractNonlinearElasticityAssembler< DIM >::GetNumNewtonIterations (  )  [inline]

Get number of Newton iterations taken in last solve.

Definition at line 871 of file AbstractNonlinearElasticityAssembler.hpp.

References AbstractNonlinearElasticityAssembler< DIM >::mNumNewtonIterations.

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

template<unsigned DIM>
void AbstractNonlinearElasticityAssembler< DIM >::PostNewtonStep ( unsigned  counter,
double  normResidual 
) [inline, protected, virtual]

This function may be overloaded by subclasses. It is called after each Newton iteration.

Parameters:
counter 
normResidual 

Definition at line 693 of file AbstractNonlinearElasticityAssembler.hpp.

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

template<unsigned DIM>
void AbstractNonlinearElasticityAssembler< DIM >::PrintLineSearchResult ( double  s,
double  residNorm 
) [inline, protected]

Print to std::cout the residual norm for this s, ie ||f(x+su)|| where f is the residual vector, x the current solution and u the update vector

Parameters:
s s
residNorm residual norm.

Definition at line 540 of file AbstractNonlinearElasticityAssembler.hpp.

Referenced by AbstractNonlinearElasticityAssembler< DIM >::UpdateSolutionUsingLineSearch().

template<unsigned DIM>
virtual std::vector<c_vector<double,DIM> >& AbstractNonlinearElasticityAssembler< DIM >::rGetDeformedPosition (  )  [protected, pure virtual]

Get the deformed position. Must be overridden in concrete derived classes.

Implemented in NonlinearElasticityAssembler< DIM >.

Referenced by AbstractNonlinearElasticityAssembler< DIM >::WriteOutput().

template<unsigned DIM>
void AbstractNonlinearElasticityAssembler< DIM >::SetFunctionalBodyForce ( c_vector< double, DIM >(*)(c_vector< double, DIM > &)  pFunction  )  [inline]

Set a function which gives body force as a function of X (undeformed position) Whatever body force was provided in the constructor will now be ignored.

Parameters:
pFunction 

Definition at line 878 of file AbstractNonlinearElasticityAssembler.hpp.

References AbstractNonlinearElasticityAssembler< DIM >::mpBodyForceFunction, and AbstractNonlinearElasticityAssembler< DIM >::mUsingBodyForceFunction.

template<unsigned DIM>
void AbstractNonlinearElasticityAssembler< DIM >::SetWriteOutput ( bool  writeOutput = true  )  [inline]

Set whether to write any output.

Parameters:
writeOutput (defaults to true)

Definition at line 886 of file AbstractNonlinearElasticityAssembler.hpp.

References AbstractNonlinearElasticityAssembler< DIM >::mOutputDirectory, and AbstractNonlinearElasticityAssembler< DIM >::mWriteOutput.

template<unsigned DIM>
void AbstractNonlinearElasticityAssembler< DIM >::Solve ( double  tol = -1.0,
unsigned  offset = 0,
unsigned  maxNumNewtonIterations = INT_MAX,
bool  quitIfNoConvergence = true 
) [inline]
template<unsigned DIM>
double AbstractNonlinearElasticityAssembler< DIM >::TakeNewtonStep (  )  [inline, protected]
template<unsigned DIM>
double AbstractNonlinearElasticityAssembler< DIM >::UpdateSolutionUsingLineSearch ( Vec  solution  )  [inline, protected]

Using the update vector (of Newton's method), choose s such that ||f(x+su)|| is most decreased, where f is the residual vector, x the current solution (mCurrentSolution) and u the update vector. This checks s=1 first (most likely to be the current solution, then 0.9, 0.8.. until ||f|| starts increasing.

Parameters:
solution The solution of the linear solve in newton's method, ie the update vector u.

Definition at line 548 of file AbstractNonlinearElasticityAssembler.hpp.

References AbstractNonlinearElasticityAssembler< DIM >::CalculateResidualNorm(), AbstractNonlinearElasticityAssembler< DIM >::ComputeResidualAndGetNorm(), AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, AbstractNonlinearElasticityAssembler< DIM >::PrintLineSearchResult(), and AbstractNonlinearElasticityAssembler< DIM >::VectorSum().

Referenced by AbstractNonlinearElasticityAssembler< DIM >::TakeNewtonStep().

template<unsigned DIM>
void AbstractNonlinearElasticityAssembler< DIM >::VectorSum ( std::vector< double > &  rX,
ReplicatableVector rY,
double  a,
std::vector< double > &  rZ 
) [inline, protected]

Simple helper function, computes Z = X + aY, where X and Z are std::vectors and Y a ReplicatableVector

Parameters:
rX X
rY Y (replicatable vector)
a a
rZ Z the returned vector

Definition at line 425 of file AbstractNonlinearElasticityAssembler.hpp.

References ReplicatableVector::GetSize().

Referenced by AbstractNonlinearElasticityAssembler< DIM >::UpdateSolutionUsingLineSearch().

template<unsigned DIM>
void AbstractNonlinearElasticityAssembler< DIM >::WriteOutput ( unsigned  counter  )  [inline]

Member Data Documentation

template<unsigned DIM>
FourthOrderTensor<DIM> AbstractNonlinearElasticityAssembler< DIM >::dTdE [protected]

Storage space for a 4th order tensor used in assembling the Jacobian (to avoid repeated memory allocation)

Definition at line 136 of file AbstractNonlinearElasticityAssembler.hpp.

Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnElement(), and AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement().

template<unsigned DIM>
const double AbstractNonlinearElasticityAssembler< DIM >::MAX_NEWTON_ABS_TOL = 1e-7 [inline, static, protected]

Maximum absolute tolerance for Newton solve.

Definition at line 59 of file AbstractNonlinearElasticityAssembler.hpp.

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

template<unsigned DIM>
c_vector<double,DIM> AbstractNonlinearElasticityAssembler< DIM >::mBodyForce [protected]
template<unsigned DIM>
std::vector<double> AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution [protected]
template<unsigned DIM>
std::vector<c_vector<double,DIM> > AbstractNonlinearElasticityAssembler< DIM >::mDeformedPosition [protected]

Deformed position: mDeformedPosition[i](j) = x_j for node i

Definition at line 142 of file AbstractNonlinearElasticityAssembler.hpp.

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

template<unsigned DIM>
double AbstractNonlinearElasticityAssembler< DIM >::mDensity [protected]

Mass density of the undeformed body (equal to the density of deformed body)

Definition at line 111 of file AbstractNonlinearElasticityAssembler.hpp.

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

template<unsigned DIM>
std::vector<c_vector<double,DIM> > AbstractNonlinearElasticityAssembler< DIM >::mFixedNodeDisplacements [protected]

The displacements of those nodes with displacement boundary conditions

Definition at line 120 of file AbstractNonlinearElasticityAssembler.hpp.

Referenced by AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), and NonlinearElasticityAssembler< DIM >::Initialise().

template<unsigned DIM>
std::vector<unsigned> AbstractNonlinearElasticityAssembler< DIM >::mFixedNodes [protected]
template<unsigned DIM>
const double AbstractNonlinearElasticityAssembler< DIM >::MIN_NEWTON_ABS_TOL = 1e-10 [inline, static, protected]

Minimum absolute tolerance for Newton solv.e

Definition at line 62 of file AbstractNonlinearElasticityAssembler.hpp.

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

template<unsigned DIM>
std::vector<AbstractIncompressibleMaterialLaw<DIM>*> AbstractNonlinearElasticityAssembler< DIM >::mMaterialLaws [protected]
template<unsigned DIM>
unsigned AbstractNonlinearElasticityAssembler< DIM >::mNumDofs [protected]

Number of degrees of freedom (eg equal to DIM*N + M if quadratic-linear bases are used, where there are N total nodes and M vertices).

Definition at line 71 of file AbstractNonlinearElasticityAssembler.hpp.

Referenced by NonlinearElasticityAssembler< DIM >::AllocateMatrixMemory(), NonlinearElasticityAssembler< DIM >::AssembleSystem(), AbstractNonlinearElasticityAssembler< DIM >::CalculateResidualNorm(), and NonlinearElasticityAssembler< DIM >::FormInitialGuess().

template<unsigned DIM>
unsigned AbstractNonlinearElasticityAssembler< DIM >::mNumNewtonIterations [protected]
template<unsigned DIM>
std::string AbstractNonlinearElasticityAssembler< DIM >::mOutputDirectory [protected]
template<unsigned DIM>
c_vector<double,DIM>(* AbstractNonlinearElasticityAssembler< DIM >::mpBodyForceFunction)(c_vector< double, DIM > &) [protected]

An optionally provided (pointer to a) function, giving body force as a function of undeformed position.

Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnElement(), and AbstractNonlinearElasticityAssembler< DIM >::SetFunctionalBodyForce().

template<unsigned DIM>
LinearSystem* AbstractNonlinearElasticityAssembler< DIM >::mpLinearSystem [protected]

The linear system which stores the matrix used for preconditioning (given the helper functions on LinearSystem it is best to use LinearSystem and use these for assembling the preconditioner, rather than just use a Mat The preconditioner is the petsc LU factorisation of

Jp = [A B] in displacement-pressure block form, [C M]

where the A, B and C are the matrices in the normal jacobian, ie

J = [A B] [C 0]

and M is the MASS MATRIX (ie integral phi_i phi_j dV, where phi_i are the pressure basis functions).

Definition at line 105 of file AbstractNonlinearElasticityAssembler.hpp.

Referenced by NonlinearElasticityAssembler< DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), NonlinearElasticityAssembler< DIM >::AssembleSystem(), AbstractNonlinearElasticityAssembler< DIM >::TakeNewtonStep(), and AbstractNonlinearElasticityAssembler< DIM >::~AbstractNonlinearElasticityAssembler().

template<unsigned DIM>
std::vector<double> AbstractNonlinearElasticityAssembler< DIM >::mPressures [protected]

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

Definition at line 148 of file AbstractNonlinearElasticityAssembler.hpp.

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

template<unsigned DIM>
c_vector<double,DIM>(* AbstractNonlinearElasticityAssembler< DIM >::mpTractionBoundaryConditionFunction)(c_vector< double, DIM > &) [protected]

An optionally provided (pointer to a) function, giving the surface traction as a function of undeformed position.

Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnBoundaryElement(), and NonlinearElasticityAssembler< DIM >::SetFunctionalTractionBoundaryCondition().

template<unsigned DIM>
std::vector<c_vector<double,DIM> > AbstractNonlinearElasticityAssembler< DIM >::mSurfaceTractions [protected]

The surface tractions (which should really be non-zero) for the boundary elements in mBoundaryElements.

Definition at line 154 of file AbstractNonlinearElasticityAssembler.hpp.

Referenced by NonlinearElasticityAssembler< DIM >::AssembleSystem(), and NonlinearElasticityAssembler< DIM >::SetSurfaceTractionBoundaryConditions().

template<unsigned DIM>
bool AbstractNonlinearElasticityAssembler< DIM >::mUsingBodyForceFunction [protected]

Whether the functional version of the body force is being used or not

Definition at line 166 of file AbstractNonlinearElasticityAssembler.hpp.

Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnElement(), and AbstractNonlinearElasticityAssembler< DIM >::SetFunctionalBodyForce().

Whether the functional version of the surface traction is being used or not

Definition at line 169 of file AbstractNonlinearElasticityAssembler.hpp.

Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnBoundaryElement(), and NonlinearElasticityAssembler< DIM >::SetFunctionalTractionBoundaryCondition().

template<unsigned DIM>
bool AbstractNonlinearElasticityAssembler< DIM >::mWriteOutput [protected]
template<unsigned DIM>
const double AbstractNonlinearElasticityAssembler< DIM >::NEWTON_REL_TOL = 1e-4 [inline, static, protected]

Relative tolerance for Newton solve.

Definition at line 65 of file AbstractNonlinearElasticityAssembler.hpp.

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


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

Generated by  doxygen 1.6.2