AbstractNonlinearElasticitySolver< DIM > Class Template Reference

#include <AbstractNonlinearElasticitySolver.hpp>

Inheritance diagram for AbstractNonlinearElasticitySolver< DIM >:

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

Collaboration graph
[legend]

List of all members.

Public Member Functions

 AbstractNonlinearElasticitySolver (QuadraticMesh< DIM > *pQuadMesh, c_vector< double, DIM > bodyForce, double density, std::string outputDirectory, std::vector< unsigned > &fixedNodes, CompressibilityType compressibilityType)
virtual ~AbstractNonlinearElasticitySolver ()
void Solve (double tol=-1.0, unsigned maxNumNewtonIterations=INT_MAX, bool quitIfNoConvergence=true)
void WriteCurrentDeformation (std::string fileName, int counterToAppend=-1)
unsigned GetNumNewtonIterations ()
void SetFunctionalBodyForce (c_vector< double, DIM >(*pFunction)(c_vector< double, DIM > &X, double t))
void SetWriteOutput (bool writeOutput=true)
void SetWriteOutputEachNewtonIteration (bool writeOutputEachNewtonIteration=true)
void SetKspAbsoluteTolerance (double kspAbsoluteTolerance)
std::vector< double > & rGetCurrentSolution ()
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 > &X, double t))
std::vector< c_vector< double,
DIM > > & 
rGetDeformedPosition ()
void SetCurrentTime (double time)
void CreateCmguiOutput ()

Protected Member Functions

virtual void AssembleSystem (bool assembleResidual, bool assembleLinearSystem)=0
void Initialise (std::vector< c_vector< double, DIM > > *pFixedNodeLocations)
void AllocateMatrixMemory ()
void ApplyBoundaryConditions (bool applyToMatrix)
void FinishAssembleSystem (bool assembleResidual, bool assembleLinearSystem)
double ComputeResidualAndGetNorm (bool allowException)
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)
virtual void ComputeStressAndStressDerivative (AbstractMaterialLaw< 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
double mKspAbsoluteTol
unsigned mNumDofs
Vec mResidualVector
Vec mLinearSystemRhsVector
Mat mJacobianMatrix
Vec mDirichletBoundaryConditionsVector
Mat mPreconditionMatrix
c_vector< double, DIM > mBodyForce
double mDensity
std::vector< unsigned > mFixedNodes
std::vector< c_vector< double,
DIM > > 
mFixedNodeDisplacements
bool mWriteOutput
std::string mOutputDirectory
OutputFileHandlermpOutputFileHandler
bool mWriteOutputEachNewtonIteration
std::vector< double > mCurrentSolution
FourthOrderTensor< DIM, DIM,
DIM, DIM > 
dTdE
unsigned mNumNewtonIterations
std::vector< c_vector< double,
DIM > > 
mDeformedPosition
std::vector< c_vector< double,
DIM > > 
mSurfaceTractions
c_vector< double, DIM >(* mpBodyForceFunction )(c_vector< double, DIM > &X, double t)
c_vector< double, DIM >(* mpTractionBoundaryConditionFunction )(c_vector< double, DIM > &X, double t)
bool mUsingBodyForceFunction
bool mUsingTractionBoundaryConditionFunction
double mCurrentTime
CompressibilityType mCompressibilityType

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 NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2
static double MAX_NEWTON_ABS_TOL = 1e-7
static double MIN_NEWTON_ABS_TOL = 1e-10
static double NEWTON_REL_TOL = 1e-4


Detailed Description

template<unsigned DIM>
class AbstractNonlinearElasticitySolver< DIM >

Abstract nonlinear elasticity solver.

Definition at line 76 of file AbstractNonlinearElasticitySolver.hpp.


Constructor & Destructor Documentation

template<unsigned DIM>
AbstractNonlinearElasticitySolver< DIM >::AbstractNonlinearElasticitySolver ( QuadraticMesh< DIM > *  pQuadMesh,
c_vector< double, DIM >  bodyForce,
double  density,
std::string  outputDirectory,
std::vector< unsigned > &  fixedNodes,
CompressibilityType  compressibilityType 
) [inline]

Constructor.

Parameters:
pQuadMesh the quadratic mesh
bodyForce Body force density (for example, acceleration due to gravity)
density density
outputDirectory output directory
fixedNodes std::vector of nodes which have a dirichlet boundary condition imposed on them
compressibilityType Should be equal to COMPRESSIBLE or INCOMPRESSIBLE (see enumeration defined at top of file) (depending on which concrete class is inheriting from this) and is only used in computing mNumDofs and allocating matrix memory.

Definition at line 1262 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mCompressibilityType, AbstractNonlinearElasticitySolver< DIM >::mNumDofs, AbstractNonlinearElasticitySolver< DIM >::mOutputDirectory, AbstractNonlinearElasticitySolver< DIM >::mpOutputFileHandler, AbstractNonlinearElasticitySolver< DIM >::mpQuadMesh, and AbstractNonlinearElasticitySolver< DIM >::mWriteOutput.

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


Member Function Documentation

template<unsigned DIM>
virtual void AbstractNonlinearElasticitySolver< DIM >::AssembleSystem ( bool  assembleResidual,
bool  assembleLinearSystem 
) [protected, pure virtual]

Assemble the residual vector and/or Jacobian matrix (using the current solution stored in mCurrentSolution, output going to mResidualVector and/or mJacobianMatrix).

Must be overridden in concrete derived classes.

Parameters:
assembleResidual A bool stating whether to assemble the residual vector.
assembleLinearSystem A bool stating whether to assemble the Jacobian matrix and the RHS vector of the linear system (which is based on the residual but could be slightly different due to the way dirichlet boundary conditions are applied to the linear system - see comments in ApplyBoundaryConditions).

Implemented in CompressibleNonlinearElasticitySolver< DIM >, and NonlinearElasticitySolver< DIM >.

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

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::Initialise ( std::vector< c_vector< double, DIM > > *  pFixedNodeLocations  )  [inline, protected]

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::AllocateMatrixMemory (  )  [inline, protected]

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::ApplyBoundaryConditions ( bool  applyToMatrix  )  [inline, protected]

Apply the Dirichlet boundary conditions to the linear system.

This will always apply the Dirichlet boundary conditions to the residual vector (basically, setting a component to the difference between the current value and the correct value).

If the boolean parameter is true, this will apply the boundary conditions to the Jacobian and the linear system RHS vector (which should be equal to the residual on entering this function). In the compressible case the boundary conditions are applied by zeroing both rows and columns of the Jacobian matrix (to maintain) symmetry, which means additional changes are needed for the RHS vector.

Parameters:
applyToMatrix Whether to apply the boundary conditions to the linear system (as well as the residual).

Definition at line 755 of file AbstractNonlinearElasticitySolver.hpp.

References PetscVecTools::AddScaledVector(), PetscMatTools::AssembleFinal(), PetscMatTools::GetMatrixRowDistributed(), AbstractNonlinearElasticitySolver< DIM >::mCompressibilityType, AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< DIM >::mDirichletBoundaryConditionsVector, AbstractNonlinearElasticitySolver< DIM >::mFixedNodeDisplacements, AbstractNonlinearElasticitySolver< DIM >::mFixedNodes, AbstractNonlinearElasticitySolver< DIM >::mJacobianMatrix, AbstractNonlinearElasticitySolver< DIM >::mLinearSystemRhsVector, AbstractNonlinearElasticitySolver< DIM >::mPreconditionMatrix, AbstractNonlinearElasticitySolver< DIM >::mResidualVector, PetscVecTools::SetElement(), PetscVecTools::Zero(), PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(), and PetscMatTools::ZeroRowsWithValueOnDiagonal().

Referenced by AbstractNonlinearElasticitySolver< DIM >::FinishAssembleSystem().

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::FinishAssembleSystem ( bool  assembleResidual,
bool  assembleLinearSystem 
) [inline, protected]

template<unsigned DIM>
double AbstractNonlinearElasticitySolver< DIM >::ComputeResidualAndGetNorm ( bool  allowException  )  [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)

Parameters:
allowException Sometimes the current solution solution will be such that the residual vector cannot be computed, as (say) the material law will throw an exception as the strains are too large. If this bool is set to true, the exception will be caught, and DBL_MAX returned as the residual norm

Definition at line 896 of file AbstractNonlinearElasticitySolver.hpp.

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

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

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

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< 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 931 of file AbstractNonlinearElasticitySolver.hpp.

References ReplicatableVector::GetSize().

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

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< 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 1102 of file AbstractNonlinearElasticitySolver.hpp.

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

template<unsigned DIM>
double AbstractNonlinearElasticitySolver< DIM >::TakeNewtonStep (  )  [inline, protected]

template<unsigned DIM>
double AbstractNonlinearElasticitySolver< 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 1110 of file AbstractNonlinearElasticitySolver.hpp.

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

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

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< 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 current newton iteration number
normResidual norm of the residual

Definition at line 1256 of file AbstractNonlinearElasticitySolver.hpp.

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

template<unsigned DIM>
virtual void AbstractNonlinearElasticitySolver< DIM >::ComputeStressAndStressDerivative ( AbstractMaterialLaw< 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 given, 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.

Parameters:
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.

Definition at line 396 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractMaterialLaw< DIM >::ComputeStressAndStressDerivative().

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

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::Solve ( double  tol = -1.0,
unsigned  maxNumNewtonIterations = INT_MAX,
bool  quitIfNoConvergence = true 
) [inline]

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::WriteCurrentDeformation ( std::string  fileName,
int  counterToAppend = -1 
) [inline]

Write the current deformation of the nodes.

Parameters:
fileName (stem)
counterToAppend append a counter to the file name
For example: WriteCurrentDeformation("solution") --> file called "solution.nodes" WriteCurrentDeformation("solution",3) --> file called "solution_3.nodes"

Definition at line 1424 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mpOutputFileHandler, AbstractNonlinearElasticitySolver< DIM >::mWriteOutput, OutputFileHandler::OpenOutputFile(), and AbstractNonlinearElasticitySolver< DIM >::rGetDeformedPosition().

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

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

Get number of Newton iterations taken in last solve.

Definition at line 1456 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mNumNewtonIterations.

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::SetFunctionalBodyForce ( c_vector< double, DIM >(*)(c_vector< double, DIM > &X, double t)  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 the function, which should be a function of space and time Note that SetCurrentTime() should be called each timestep if the force changes with time

Definition at line 1463 of file AbstractNonlinearElasticitySolver.hpp.

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

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

Set whether to write any output.

Parameters:
writeOutput (defaults to true)

Definition at line 1471 of file AbstractNonlinearElasticitySolver.hpp.

References EXCEPTION, AbstractNonlinearElasticitySolver< DIM >::mOutputDirectory, and AbstractNonlinearElasticitySolver< DIM >::mWriteOutput.

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::SetWriteOutputEachNewtonIteration ( bool  writeOutputEachNewtonIteration = true  )  [inline]

By default only the original and converged solutions are written. Call this by get node positions written after every Newton step (mostly for debugging).

Parameters:
writeOutputEachNewtonIteration whether to write each iteration

Definition at line 487 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mWriteOutputEachNewtonIteration.

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::SetKspAbsoluteTolerance ( double  kspAbsoluteTolerance  )  [inline]

Set the absolute tolerance to be used when solving the linear system. If this is not called a relative tolerance is used.

Parameters:
kspAbsoluteTolerance the tolerance

Definition at line 497 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mKspAbsoluteTol.

template<unsigned DIM>
std::vector<double>& AbstractNonlinearElasticitySolver< DIM >::rGetCurrentSolution (  )  [inline]

Get the current solution vector (advanced use only - for getting the deformed position use rGetDeformedPosition()

Definition at line 507 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution.

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< 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.

Parameters:
rBoundaryElements the boundary elements
rSurfaceTractions the corresponding tractions

Definition at line 1481 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mBoundaryElements, and AbstractNonlinearElasticitySolver< DIM >::mSurfaceTractions.

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::SetFunctionalTractionBoundaryCondition ( std::vector< BoundaryElement< DIM-1, DIM > * >  rBoundaryElements,
c_vector< double, DIM >(*)(c_vector< double, DIM > &X, double t)  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.

Parameters:
rBoundaryElements the boundary elements
pFunction the function, which should be a function of space and time. Note that SetCurrentTime() should be called each timestep if the traction changes with time

Definition at line 1491 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mBoundaryElements, AbstractNonlinearElasticitySolver< DIM >::mpTractionBoundaryConditionFunction, and AbstractNonlinearElasticitySolver< DIM >::mUsingTractionBoundaryConditionFunction.

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

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::SetCurrentTime ( double  time  )  [inline]

This solver is for static problems, however the body force or surface tractions could be a function of time. This method is for setting the time.

Parameters:
time current time

Definition at line 543 of file AbstractNonlinearElasticitySolver.hpp.

References AbstractNonlinearElasticitySolver< DIM >::mCurrentTime.

template<unsigned DIM>
void AbstractNonlinearElasticitySolver< DIM >::CreateCmguiOutput (  )  [inline]


Member Data Documentation

template<unsigned DIM>
const size_t AbstractNonlinearElasticitySolver< DIM >::NUM_VERTICES_PER_ELEMENT = DIM+1 [static, protected]

Number of vertices per element

Reimplemented in CompressibleNonlinearElasticitySolver< DIM >, and NonlinearElasticitySolver< DIM >.

Definition at line 80 of file AbstractNonlinearElasticitySolver.hpp.

template<unsigned DIM>
const size_t AbstractNonlinearElasticitySolver< DIM >::NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2 [static, protected]

Number of nodes per element

Reimplemented in CompressibleNonlinearElasticitySolver< DIM >, and NonlinearElasticitySolver< DIM >.

Definition at line 82 of file AbstractNonlinearElasticitySolver.hpp.

template<unsigned DIM>
const size_t AbstractNonlinearElasticitySolver< DIM >::NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2 [static, protected]

Number of nodes per boundary element

Reimplemented in CompressibleNonlinearElasticitySolver< DIM >, and NonlinearElasticitySolver< DIM >.

Definition at line 84 of file AbstractNonlinearElasticitySolver.hpp.

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

Maximum absolute tolerance for Newton solve. The Newton solver uses the absolute tolerance corresponding to the specified relative tolerance, but has a max and min allowable absolute tolerance. (ie: if max_abs = 1e-7, min_abs = 1e-10, rel=1e-4: then if the norm of the initial_residual (=a) is 1e-2, it will solve with tolerance 1e-7; if a=1e-5, it will solve with tolerance 1e-9; a=1e-9, it will solve with tolerance 1e-10.

Definition at line 92 of file AbstractNonlinearElasticitySolver.hpp.

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

template<unsigned DIM>
double AbstractNonlinearElasticitySolver< DIM >::MIN_NEWTON_ABS_TOL = 1e-10 [inline, static, protected]

Minimum absolute tolerance for Newton solve. See documentation for MAX_NEWTON_ABS_TOL.

Definition at line 95 of file AbstractNonlinearElasticitySolver.hpp.

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

template<unsigned DIM>
double AbstractNonlinearElasticitySolver< DIM >::NEWTON_REL_TOL = 1e-4 [inline, static, protected]

Relative tolerance for Newton solve. See documentation for MAX_NEWTON_ABS_TOL.

Definition at line 98 of file AbstractNonlinearElasticitySolver.hpp.

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

template<unsigned DIM>
QuadraticMesh<DIM>* AbstractNonlinearElasticitySolver< DIM >::mpQuadMesh [protected]

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

template<unsigned DIM>
GaussianQuadratureRule<DIM>* AbstractNonlinearElasticitySolver< DIM >::mpQuadratureRule [protected]

template<unsigned DIM>
GaussianQuadratureRule<DIM-1>* AbstractNonlinearElasticitySolver< DIM >::mpBoundaryQuadratureRule [protected]

template<unsigned DIM>
double AbstractNonlinearElasticitySolver< DIM >::mKspAbsoluteTol [protected]

Absolute tolerance for linear systems. Can be set by calling SetKspAbsoluteTolerances(), but default to -1, in which case a relative tolerance is used.

Definition at line 120 of file AbstractNonlinearElasticitySolver.hpp.

Referenced by AbstractNonlinearElasticitySolver< DIM >::SetKspAbsoluteTolerance(), and AbstractNonlinearElasticitySolver< DIM >::TakeNewtonStep().

template<unsigned DIM>
unsigned AbstractNonlinearElasticitySolver< DIM >::mNumDofs [protected]

template<unsigned DIM>
Vec AbstractNonlinearElasticitySolver< DIM >::mResidualVector [protected]

template<unsigned DIM>
Vec AbstractNonlinearElasticitySolver< DIM >::mLinearSystemRhsVector [protected]

The RHS side in the linear system that is solved each Newton iteration. Since Newton's method is Ju = f, where J is the Jacobian, u the (negative of the) update and f the residual, it might seem necessary to store this as well as the residual. However, when applying Dirichlet boundary conditions in the compressible case, we alter the rows of the matrix, but also alter the columns in order to maintain symmetry. This requires making further changes to the right-hand vector, meaning that it no longer properly represents the residual. Hence, we have to use two vectors.

Overall, this can be represents as

  • compute residual f
  • compute Jacobian J
  • apply BCs to f.
  • alter the linear system from Ju=f to (J*)u=f* which enforces the dirichlet boundary conditions but enforces them symmetrically.

mLinearSystemRhsVector represents f*.

Definition at line 152 of file AbstractNonlinearElasticitySolver.hpp.

Referenced by AbstractNonlinearElasticitySolver< DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticitySolver< DIM >::ApplyBoundaryConditions(), AbstractNonlinearElasticitySolver< DIM >::FinishAssembleSystem(), AbstractNonlinearElasticitySolver< DIM >::TakeNewtonStep(), and AbstractNonlinearElasticitySolver< DIM >::~AbstractNonlinearElasticitySolver().

template<unsigned DIM>
Mat AbstractNonlinearElasticitySolver< DIM >::mJacobianMatrix [protected]

template<unsigned DIM>
Vec AbstractNonlinearElasticitySolver< DIM >::mDirichletBoundaryConditionsVector [protected]

template<unsigned DIM>
Mat AbstractNonlinearElasticitySolver< DIM >::mPreconditionMatrix [protected]

Precondition matrix for the linear system

In the incompressible case: 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 183 of file AbstractNonlinearElasticitySolver.hpp.

Referenced by AbstractNonlinearElasticitySolver< DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticitySolver< DIM >::ApplyBoundaryConditions(), NonlinearElasticitySolver< DIM >::AssembleSystem(), CompressibleNonlinearElasticitySolver< DIM >::AssembleSystem(), AbstractNonlinearElasticitySolver< DIM >::FinishAssembleSystem(), AbstractNonlinearElasticitySolver< DIM >::TakeNewtonStep(), and AbstractNonlinearElasticitySolver< DIM >::~AbstractNonlinearElasticitySolver().

template<unsigned DIM>
c_vector<double,DIM> AbstractNonlinearElasticitySolver< DIM >::mBodyForce [protected]

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

Mass density of the undeformed body (equal to the density of deformed body in the incompressible case)

Definition at line 189 of file AbstractNonlinearElasticitySolver.hpp.

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

template<unsigned DIM>
std::vector<unsigned> AbstractNonlinearElasticitySolver< DIM >::mFixedNodes [protected]

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

The displacements of those nodes with displacement boundary conditions

Definition at line 195 of file AbstractNonlinearElasticitySolver.hpp.

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

template<unsigned DIM>
bool AbstractNonlinearElasticitySolver< DIM >::mWriteOutput [protected]

template<unsigned DIM>
std::string AbstractNonlinearElasticitySolver< DIM >::mOutputDirectory [protected]

template<unsigned DIM>
OutputFileHandler* AbstractNonlinearElasticitySolver< DIM >::mpOutputFileHandler [protected]

template<unsigned DIM>
bool AbstractNonlinearElasticitySolver< DIM >::mWriteOutputEachNewtonIteration [protected]

By default only the initial and final solutions are written. However, we may want to write the solutions after every Newton iteration, in which case the following should be set to true

Definition at line 209 of file AbstractNonlinearElasticitySolver.hpp.

Referenced by AbstractNonlinearElasticitySolver< DIM >::SetWriteOutputEachNewtonIteration(), and AbstractNonlinearElasticitySolver< DIM >::Solve().

template<unsigned DIM>
std::vector<double> AbstractNonlinearElasticitySolver< DIM >::mCurrentSolution [protected]

template<unsigned DIM>
FourthOrderTensor<DIM,DIM,DIM,DIM> AbstractNonlinearElasticitySolver< DIM >::dTdE [protected]

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

Definition at line 223 of file AbstractNonlinearElasticitySolver.hpp.

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

template<unsigned DIM>
unsigned AbstractNonlinearElasticitySolver< DIM >::mNumNewtonIterations [protected]

template<unsigned DIM>
std::vector<c_vector<double,DIM> > AbstractNonlinearElasticitySolver< DIM >::mDeformedPosition [protected]

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

Definition at line 229 of file AbstractNonlinearElasticitySolver.hpp.

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

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

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

Definition at line 235 of file AbstractNonlinearElasticitySolver.hpp.

Referenced by NonlinearElasticitySolver< DIM >::AssembleSystem(), CompressibleNonlinearElasticitySolver< DIM >::AssembleSystem(), and AbstractNonlinearElasticitySolver< DIM >::SetSurfaceTractionBoundaryConditions().

template<unsigned DIM>
c_vector<double,DIM>(* AbstractNonlinearElasticitySolver< DIM >::mpBodyForceFunction)(c_vector< double, DIM > &X, double t) [protected]

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

Referenced by NonlinearElasticitySolver< DIM >::AssembleOnElement(), CompressibleNonlinearElasticitySolver< DIM >::AssembleOnElement(), and AbstractNonlinearElasticitySolver< DIM >::SetFunctionalBodyForce().

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

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

template<unsigned DIM>
bool AbstractNonlinearElasticitySolver< DIM >::mUsingTractionBoundaryConditionFunction [protected]

template<unsigned DIM>
double AbstractNonlinearElasticitySolver< DIM >::mCurrentTime [protected]

template<unsigned DIM>
CompressibilityType AbstractNonlinearElasticitySolver< DIM >::mCompressibilityType [protected]


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

Generated on Tue May 31 14:32:21 2011 for Chaste by  doxygen 1.5.5