#include <AbstractNonlinearElasticityAssembler.hpp>
Inherited by NonlinearElasticityAssembler< DIM >.
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 |
LinearSystem * | mpLinearSystem |
LinearSystem * | mpPreconditionMatrixLinearSystem |
c_vector< double, DIM > | mBodyForce |
double | mDensity |
std::string | mOutputDirectory |
std::vector< unsigned > | mFixedNodes |
std::vector< c_vector< double, DIM > > | mFixedNodeDisplacements |
bool | mWriteOutput |
std::vector< double > | mCurrentSolution |
FourthOrderTensor< DIM > | dTdE |
unsigned | mNumNewtonIterations |
std::vector< c_vector< double, DIM > > | mDeformedPosition |
std::vector< double > | mPressures |
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 |
Abstract nonlinear elasticity assembler.
Definition at line 54 of file AbstractNonlinearElasticityAssembler.hpp.
AbstractNonlinearElasticityAssembler< DIM >::AbstractNonlinearElasticityAssembler | ( | unsigned | numDofs, | |
AbstractIncompressibleMaterialLaw< DIM > * | pMaterialLaw, | |||
c_vector< double, DIM > | bodyForce, | |||
double | density, | |||
std::string | outputDirectory, | |||
std::vector< unsigned > & | fixedNodes | |||
) | [inline] |
Constructor.
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.
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.
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.
AbstractNonlinearElasticityAssembler< DIM >::~AbstractNonlinearElasticityAssembler | ( | ) | [inline, virtual] |
Destructor.
Definition at line 754 of file AbstractNonlinearElasticityAssembler.hpp.
References AbstractNonlinearElasticityAssembler< DIM >::mpLinearSystem, and AbstractNonlinearElasticityAssembler< DIM >::mpPreconditionMatrixLinearSystem.
void AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions | ( | bool | applyToMatrix | ) | [inline, protected] |
Apply the Dirichlet boundary conditions to the linear system.
applyToMatrix |
Definition at line 347 of file AbstractNonlinearElasticityAssembler.hpp.
References AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution, AbstractNonlinearElasticityAssembler< DIM >::mFixedNodeDisplacements, AbstractNonlinearElasticityAssembler< DIM >::mFixedNodes, AbstractNonlinearElasticityAssembler< DIM >::mpLinearSystem, AbstractNonlinearElasticityAssembler< DIM >::mpPreconditionMatrixLinearSystem, LinearSystem::SetRhsVectorElement(), and LinearSystem::ZeroMatrixRowsWithValueOnDiagonal().
Referenced by NonlinearElasticityAssembler< DIM >::AssembleSystem().
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.
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().
double AbstractNonlinearElasticityAssembler< DIM >::CalculateResidualNorm | ( | ) | [inline, protected] |
Calculate |r|_2 / length(r), where r is the current residual vector
Definition at line 417 of file AbstractNonlinearElasticityAssembler.hpp.
References AbstractNonlinearElasticityAssembler< DIM >::mNumDofs, AbstractNonlinearElasticityAssembler< DIM >::mpLinearSystem, and LinearSystem::rGetRhsVector().
Referenced by AbstractNonlinearElasticityAssembler< DIM >::ComputeResidualAndGetNorm(), and AbstractNonlinearElasticityAssembler< DIM >::UpdateSolutionUsingLineSearch().
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().
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 >.
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().
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.
counter | ||
normResidual |
Definition at line 693 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by AbstractNonlinearElasticityAssembler< DIM >::Solve().
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
s | s | |
residNorm | residual norm. |
Definition at line 540 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by AbstractNonlinearElasticityAssembler< DIM >::UpdateSolutionUsingLineSearch().
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().
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.
pFunction |
Definition at line 878 of file AbstractNonlinearElasticityAssembler.hpp.
References AbstractNonlinearElasticityAssembler< DIM >::mpBodyForceFunction, and AbstractNonlinearElasticityAssembler< DIM >::mUsingBodyForceFunction.
void AbstractNonlinearElasticityAssembler< DIM >::SetWriteOutput | ( | bool | writeOutput = true |
) | [inline] |
Set whether to write any output.
writeOutput | (defaults to true) |
Definition at line 886 of file AbstractNonlinearElasticityAssembler.hpp.
References AbstractNonlinearElasticityAssembler< DIM >::mOutputDirectory, and AbstractNonlinearElasticityAssembler< DIM >::mWriteOutput.
void AbstractNonlinearElasticityAssembler< DIM >::Solve | ( | double | tol = -1.0 , |
|
unsigned | offset = 0 , |
|||
unsigned | maxNumNewtonIterations = INT_MAX , |
|||
bool | quitIfNoConvergence = true | |||
) | [inline] |
Solve the problem.
tol | (defaults to -1.0) | |
offset | (defaults to 0) | |
maxNumNewtonIterations | (defaults to INT_MAX) | |
quitIfNoConvergence | (defaults to true) |
Definition at line 762 of file AbstractNonlinearElasticityAssembler.hpp.
References AbstractNonlinearElasticityAssembler< DIM >::ComputeResidualAndGetNorm(), AbstractNonlinearElasticityAssembler< DIM >::MAX_NEWTON_ABS_TOL, AbstractNonlinearElasticityAssembler< DIM >::MIN_NEWTON_ABS_TOL, AbstractNonlinearElasticityAssembler< DIM >::mNumNewtonIterations, AbstractNonlinearElasticityAssembler< DIM >::mWriteOutput, AbstractNonlinearElasticityAssembler< DIM >::NEWTON_REL_TOL, AbstractNonlinearElasticityAssembler< DIM >::PostNewtonStep(), AbstractNonlinearElasticityAssembler< DIM >::TakeNewtonStep(), and AbstractNonlinearElasticityAssembler< DIM >::WriteOutput().
double AbstractNonlinearElasticityAssembler< DIM >::TakeNewtonStep | ( | ) | [inline, protected] |
Take one newton step, by solving the linear system -Ju=f, (J the jacobian, f the residual, u the update), and picking s such that a_new = a_old + su (a the current solution) such |f(a)| is the smallest.
Definition at line 440 of file AbstractNonlinearElasticityAssembler.hpp.
References AbstractNonlinearElasticityAssembler< DIM >::AssembleSystem(), GenericEventHandler< 7, MechanicsEventHandler >::BeginEvent(), GenericEventHandler< 7, MechanicsEventHandler >::EndEvent(), PetscTools::GetMyRank(), AbstractNonlinearElasticityAssembler< DIM >::mpLinearSystem, AbstractNonlinearElasticityAssembler< DIM >::mpPreconditionMatrixLinearSystem, Timer::PrintAndReset(), Timer::Reset(), LinearSystem::rGetLhsMatrix(), LinearSystem::rGetRhsVector(), and AbstractNonlinearElasticityAssembler< DIM >::UpdateSolutionUsingLineSearch().
Referenced by AbstractNonlinearElasticityAssembler< DIM >::Solve().
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.
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().
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
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().
void AbstractNonlinearElasticityAssembler< DIM >::WriteOutput | ( | unsigned | counter | ) | [inline] |
Write the current solution for the file outputdir/solution_[counter].nodes
counter |
Definition at line 843 of file AbstractNonlinearElasticityAssembler.hpp.
References AbstractNonlinearElasticityAssembler< DIM >::mOutputDirectory, AbstractNonlinearElasticityAssembler< DIM >::mWriteOutput, OutputFileHandler::OpenOutputFile(), and AbstractNonlinearElasticityAssembler< DIM >::rGetDeformedPosition().
Referenced by AbstractNonlinearElasticityAssembler< DIM >::Solve().
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().
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().
c_vector<double,DIM> AbstractNonlinearElasticityAssembler< DIM >::mBodyForce [protected] |
Body force vector
Definition at line 108 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AssembleOnElement().
std::vector<double> AbstractNonlinearElasticityAssembler< DIM >::mCurrentSolution [protected] |
The current solution, in the form (in 2d) [u1 v1 u2 v2 ... uN vN p1 p2 .. pM] where there are N total nodes and M vertices
Definition at line 130 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), NonlinearElasticityAssembler< DIM >::AssembleOnElement(), AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), NonlinearElasticityAssembler< DIM >::AssembleSystem(), NonlinearElasticityAssembler< DIM >::FormInitialGuess(), NonlinearElasticityAssembler< DIM >::rGetDeformedPosition(), NonlinearElasticityAssembler< DIM >::rGetPressures(), and AbstractNonlinearElasticityAssembler< DIM >::UpdateSolutionUsingLineSearch().
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().
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().
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().
std::vector<unsigned> AbstractNonlinearElasticityAssembler< DIM >::mFixedNodes [protected] |
All nodes (including non-vertices) which are fixed
Definition at line 117 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), and NonlinearElasticityAssembler< DIM >::Initialise().
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().
std::vector<AbstractIncompressibleMaterialLaw<DIM>*> AbstractNonlinearElasticityAssembler< 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 78 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by AbstractNonlinearElasticityAssembler< DIM >::AbstractNonlinearElasticityAssembler(), NonlinearElasticityAssembler< DIM >::AssembleOnElement(), AbstractCardiacMechanicsAssembler< DIM >::AssembleOnElement(), NonlinearElasticityAssembler< DIM >::FormInitialGuess(), and AbstractCardiacMechanicsAssembler< DIM >::~AbstractCardiacMechanicsAssembler().
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().
unsigned AbstractNonlinearElasticityAssembler< DIM >::mNumNewtonIterations [protected] |
Number of Newton iterations taken in last solve
Definition at line 139 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by AbstractNonlinearElasticityAssembler< DIM >::GetNumNewtonIterations(), and AbstractNonlinearElasticityAssembler< DIM >::Solve().
std::string AbstractNonlinearElasticityAssembler< DIM >::mOutputDirectory [protected] |
Where to write output, relative to CHASTE_TESTOUTPUT
Definition at line 114 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by AbstractNonlinearElasticityAssembler< DIM >::AbstractNonlinearElasticityAssembler(), AbstractNonlinearElasticityAssembler< DIM >::SetWriteOutput(), and AbstractNonlinearElasticityAssembler< DIM >::WriteOutput().
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().
LinearSystem* AbstractNonlinearElasticityAssembler< DIM >::mpLinearSystem [protected] |
The linear system where we store all residual vectors which are calculated and the Jacobian. Note we don't actually call Solve but solve using Petsc methods explicitly (in order to easily set number of restarts etc).
Definition at line 85 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by NonlinearElasticityAssembler< DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), NonlinearElasticityAssembler< DIM >::AssembleSystem(), AbstractNonlinearElasticityAssembler< DIM >::CalculateResidualNorm(), AbstractNonlinearElasticityAssembler< DIM >::TakeNewtonStep(), and AbstractNonlinearElasticityAssembler< DIM >::~AbstractNonlinearElasticityAssembler().
LinearSystem* AbstractNonlinearElasticityAssembler< DIM >::mpPreconditionMatrixLinearSystem [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().
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().
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().
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().
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().
bool AbstractNonlinearElasticityAssembler< DIM >::mUsingTractionBoundaryConditionFunction [protected] |
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().
bool AbstractNonlinearElasticityAssembler< DIM >::mWriteOutput [protected] |
Whether to write any output
Definition at line 123 of file AbstractNonlinearElasticityAssembler.hpp.
Referenced by AbstractNonlinearElasticityAssembler< DIM >::AbstractNonlinearElasticityAssembler(), AbstractNonlinearElasticityAssembler< DIM >::SetWriteOutput(), AbstractNonlinearElasticityAssembler< DIM >::Solve(), and AbstractNonlinearElasticityAssembler< DIM >::WriteOutput().
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().