#include <AbstractNonlinearElasticitySolver.hpp>
Public Member Functions | |
AbstractNonlinearElasticitySolver (QuadraticMesh< DIM > *pQuadMesh, c_vector< double, DIM > bodyForce, double density, std::string outputDirectory, std::vector< unsigned > &fixedNodes) | |
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 assembleJacobian)=0 |
void | Initialise (std::vector< c_vector< double, DIM > > *pFixedNodeLocations) |
void | AllocateMatrixMemory () |
void | ApplyBoundaryConditions (bool applyToMatrix) |
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) |
Protected Attributes | |
QuadraticMesh< DIM > * | mpQuadMesh |
std::vector< BoundaryElement < DIM-1, DIM > * > | mBoundaryElements |
GaussianQuadratureRule< DIM > * | mpQuadratureRule |
GaussianQuadratureRule< DIM-1 > * | mpBoundaryQuadratureRule |
double | mKspAbsoluteTol |
unsigned | mNumDofs |
LinearSystem * | mpLinearSystem |
LinearSystem * | mpPreconditionMatrixLinearSystem |
c_vector< double, DIM > | mBodyForce |
double | mDensity |
std::vector< unsigned > | mFixedNodes |
std::vector< c_vector< double, DIM > > | mFixedNodeDisplacements |
bool | mWriteOutput |
std::string | mOutputDirectory |
OutputFileHandler * | mpOutputFileHandler |
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 |
Static Protected Attributes | |
static double | MAX_NEWTON_ABS_TOL = 1e-7 |
static double | MIN_NEWTON_ABS_TOL = 1e-10 |
static double | NEWTON_REL_TOL = 1e-4 |
Definition at line 73 of file AbstractNonlinearElasticitySolver.hpp.
AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AbstractNonlinearElasticitySolver | ( | QuadraticMesh< DIM > * | pQuadMesh, | |
c_vector< double, DIM > | bodyForce, | |||
double | density, | |||
std::string | outputDirectory, | |||
std::vector< unsigned > & | fixedNodes | |||
) | [inline] |
Constructor.
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 |
Definition at line 1001 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mNumDofs, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mOutputDirectory, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpOutputFileHandler, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadMesh, and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mWriteOutput.
AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::~AbstractNonlinearElasticitySolver | ( | ) | [inline, virtual] |
Destructor.
Definition at line 1045 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpBoundaryQuadratureRule, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpLinearSystem, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpOutputFileHandler, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpPreconditionMatrixLinearSystem, and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadratureRule.
virtual void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 CompressibleNonlinearElasticitySolver< DIM >, and NonlinearElasticitySolver< DIM >.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ComputeResidualAndGetNorm(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::TakeNewtonStep().
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise | ( | std::vector< c_vector< double, DIM > > * | pFixedNodeLocations | ) | [inline, protected] |
Initialise the solver.
pFixedNodeLocations |
Definition at line 465 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mFixedNodeDisplacements, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mFixedNodes, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mNumDofs, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpBoundaryQuadratureRule, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadMesh, and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadratureRule.
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AllocateMatrixMemory | ( | ) | [inline, protected] |
Allocates memory for the Jacobian and preconditioner matrices (larger number of non-zeros per row than with say linear problems)
Definition at line 507 of file AbstractNonlinearElasticitySolver.hpp.
References LinearSystem::GetOwnershipRange(), PetscTools::IsSequential(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mNumDofs, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpLinearSystem, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpPreconditionMatrixLinearSystem, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadMesh, and LinearSystem::rGetLhsMatrix().
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise().
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ApplyBoundaryConditions | ( | bool | applyToMatrix | ) | [inline, protected] |
Apply the Dirichlet boundary conditions to the linear system.
applyToMatrix | Whether to apply the boundary conditions to the matrix (as well as the vector) |
Definition at line 606 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mFixedNodeDisplacements, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mFixedNodes, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpLinearSystem, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpPreconditionMatrixLinearSystem, LinearSystem::SetRhsVectorElement(), and LinearSystem::ZeroMatrixRowsWithValueOnDiagonal().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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)
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 647 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AssembleSystem(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CalculateResidualNorm().
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::UpdateSolutionUsingLineSearch().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CalculateResidualNorm | ( | ) | [inline, protected] |
Calculate |r|_2 / length(r), where r is the current residual vector
Definition at line 674 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mNumDofs, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpLinearSystem, and LinearSystem::rGetRhsVector().
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ComputeResidualAndGetNorm(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::UpdateSolutionUsingLineSearch().
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 682 of file AbstractNonlinearElasticitySolver.hpp.
References ReplicatableVector::GetSize().
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::UpdateSolutionUsingLineSearch().
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 841 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::UpdateSolutionUsingLineSearch().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 697 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AssembleSystem(), GenericEventHandler< 7, MechanicsEventHandler >::BeginEvent(), GenericEventHandler< 7, MechanicsEventHandler >::EndEvent(), EXCEPTION, PetscTools::GetMyRank(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mKspAbsoluteTol, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpLinearSystem, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpPreconditionMatrixLinearSystem, Timer::PrintAndReset(), Timer::Reset(), LinearSystem::rGetLhsMatrix(), LinearSystem::rGetRhsVector(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::UpdateSolutionUsingLineSearch().
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 849 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CalculateResidualNorm(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ComputeResidualAndGetNorm(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::PrintLineSearchResult(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::VectorSum().
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::TakeNewtonStep().
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::PostNewtonStep | ( | unsigned | counter, | |
double | normResidual | |||
) | [inline, protected, virtual] |
This function may be overloaded by subclasses. It is called after each Newton iteration.
counter | current newton iteration number | |
normResidual | norm of the residual |
Definition at line 995 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve().
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve | ( | double | tol = -1.0 , |
|
unsigned | maxNumNewtonIterations = INT_MAX , |
|||
bool | quitIfNoConvergence = true | |||
) | [inline] |
Solve the problem.
tol | tolerance used in Newton solve (defaults to -1.0) | |
maxNumNewtonIterations | (defaults to INT_MAX) | |
quitIfNoConvergence | (defaults to true) |
Definition at line 1059 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ComputeResidualAndGetNorm(), EXCEPTION, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::MAX_NEWTON_ABS_TOL, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::MIN_NEWTON_ABS_TOL, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mNumNewtonIterations, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mWriteOutputEachNewtonIteration, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::NEWTON_REL_TOL, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::PostNewtonStep(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::TakeNewtonStep(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::WriteCurrentDeformation().
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::WriteCurrentDeformation | ( | std::string | fileName, | |
int | counterToAppend = -1 | |||
) | [inline] |
Write the current deformation of the nodes.
fileName | (stem) | |
counterToAppend | append a counter to the file name |
Definition at line 1145 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpOutputFileHandler, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mWriteOutput, OutputFileHandler::OpenOutputFile(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::rGetDeformedPosition().
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve().
unsigned AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::GetNumNewtonIterations | ( | ) | [inline] |
Get number of Newton iterations taken in last solve.
Definition at line 1177 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mNumNewtonIterations.
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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.
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 1184 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpBodyForceFunction, and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mUsingBodyForceFunction.
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetWriteOutput | ( | bool | writeOutput = true |
) | [inline] |
Set whether to write any output.
writeOutput | (defaults to true) |
Definition at line 1192 of file AbstractNonlinearElasticitySolver.hpp.
References EXCEPTION, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mOutputDirectory, and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mWriteOutput.
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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).
writeOutputEachNewtonIteration | whether to write each iteration |
Definition at line 391 of file AbstractNonlinearElasticitySolver.hpp.
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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.
kspAbsoluteTolerance | the tolerance |
Definition at line 401 of file AbstractNonlinearElasticitySolver.hpp.
std::vector<double>& AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::rGetCurrentSolution | ( | ) | [inline] |
Get the current solution vector (advanced use only - for getting the deformed position use rGetDeformedPosition()
Definition at line 411 of file AbstractNonlinearElasticitySolver.hpp.
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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.
rBoundaryElements | the boundary elements | |
rSurfaceTractions | the corresponding tractions |
Definition at line 1202 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mBoundaryElements, and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mSurfaceTractions.
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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.
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 1212 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mBoundaryElements, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpTractionBoundaryConditionFunction, and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mUsingTractionBoundaryConditionFunction.
std::vector< c_vector< double, DIM > > & AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::rGetDeformedPosition | ( | ) | [inline] |
Get the deformed position. Note returnvalue[i](j) = x_j for node i.
Definition at line 1223 of file AbstractNonlinearElasticitySolver.hpp.
References AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mCurrentSolution, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mDeformedPosition, and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadMesh.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CreateCmguiOutput(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::WriteCurrentDeformation().
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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.
time | current time |
Definition at line 447 of file AbstractNonlinearElasticitySolver.hpp.
void AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CreateCmguiOutput | ( | ) | [inline] |
Convert the output to Cmgui format (placed in a folder called cmgui in the output directory). Writes the original mesh as solution_0.exnode and the (current) solution as solution_1.exnode
Definition at line 1238 of file AbstractNonlinearElasticitySolver.hpp.
References EXCEPTION, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mOutputDirectory, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadMesh, AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::rGetDeformedPosition(), CmguiDeformedSolutionsWriter< DIM >::WriteCmguiScript(), CmguiDeformedSolutionsWriter< DIM >::WriteDeformationPositions(), and CmguiDeformedSolutionsWriter< DIM >::WriteInitialMesh().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 82 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 85 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::NEWTON_REL_TOL = 1e-4 [inline, static, protected] |
Relative tolerance for Newton solve. See documentation for MAX_NEWTON_ABS_TOL.
Definition at line 88 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve().
QuadraticMesh<DIM>* AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadMesh [protected] |
The mesh to be solved on. Requires 6 nodes per triangle (or 10 per tetrahedron) as quadratic bases are used.
Definition at line 95 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AbstractNonlinearElasticitySolver(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CreateCmguiOutput(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::rGetDeformedPosition().
std::vector<BoundaryElement<DIM-1,DIM>*> AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mBoundaryElements [protected] |
Boundary elements with (non-zero) surface tractions defined on them
Definition at line 98 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetFunctionalTractionBoundaryCondition(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetSurfaceTractionBoundaryConditions().
GaussianQuadratureRule<DIM>* AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpQuadratureRule [protected] |
Gaussian quadrature rule
Definition at line 101 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::~AbstractNonlinearElasticitySolver().
GaussianQuadratureRule<DIM-1>* AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpBoundaryQuadratureRule [protected] |
Boundary Gaussian quadrature rule
Definition at line 104 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::~AbstractNonlinearElasticitySolver().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 110 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBLE, DIM >::SetKspAbsoluteTolerance(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::TakeNewtonStep().
unsigned AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mNumDofs [protected] |
Number of degrees of freedom (equal to, in the incompressible case: DIM*N + M if quadratic-linear bases are used, where there are N total nodes and M vertices; or DIM*N in the compressible case).
Definition at line 117 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AbstractNonlinearElasticitySolver(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CalculateResidualNorm(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise().
LinearSystem* AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 124 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ApplyBoundaryConditions(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CalculateResidualNorm(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::TakeNewtonStep(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::~AbstractNonlinearElasticitySolver().
LinearSystem* AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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.
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 147 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AllocateMatrixMemory(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ApplyBoundaryConditions(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::TakeNewtonStep(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::~AbstractNonlinearElasticitySolver().
c_vector<double,DIM> AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mBodyForce [protected] |
Body force vector
Definition at line 150 of file AbstractNonlinearElasticitySolver.hpp.
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mDensity [protected] |
Mass density of the undeformed body (equal to the density of deformed body in the incompressible case)
Definition at line 153 of file AbstractNonlinearElasticitySolver.hpp.
std::vector<unsigned> AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mFixedNodes [protected] |
All nodes (including non-vertices) which are fixed
Definition at line 156 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ApplyBoundaryConditions(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise().
std::vector<c_vector<double,DIM> > AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mFixedNodeDisplacements [protected] |
The displacements of those nodes with displacement boundary conditions
Definition at line 159 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ApplyBoundaryConditions(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise().
bool AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mWriteOutput [protected] |
Whether to write any output
Definition at line 162 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AbstractNonlinearElasticitySolver(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetWriteOutput(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::WriteCurrentDeformation().
std::string AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mOutputDirectory [protected] |
Where to write output, relative to CHASTE_TESTOUTPUT
Definition at line 165 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AbstractNonlinearElasticitySolver(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::CreateCmguiOutput(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetWriteOutput().
OutputFileHandler* AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpOutputFileHandler [protected] |
Output file handler
Definition at line 168 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::AbstractNonlinearElasticitySolver(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::WriteCurrentDeformation(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::~AbstractNonlinearElasticitySolver().
bool AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 173 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBLE, DIM >::SetWriteOutputEachNewtonIteration(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve().
std::vector<double> AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mCurrentSolution [protected] |
The current solution, in the form (assuming 2d): Incompressible problem: [u1 v1 u2 v2 ... uN vN p1 p2 .. pM] Compressible problem: [u1 v1 u2 v2 ... uN vN] where there are N total nodes and M vertices
Definition at line 182 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::ApplyBoundaryConditions(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Initialise(), AbstractNonlinearElasticitySolver< COMPRESSIBLE, DIM >::rGetCurrentSolution(), AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::rGetDeformedPosition(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::UpdateSolutionUsingLineSearch().
FourthOrderTensor<DIM,DIM,DIM,DIM> AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::dTdE [protected] |
Storage space for a 4th order tensor used in assembling the Jacobian (to avoid repeated memory allocation)
Definition at line 187 of file AbstractNonlinearElasticitySolver.hpp.
unsigned AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mNumNewtonIterations [protected] |
Number of Newton iterations taken in last solve
Definition at line 190 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::GetNumNewtonIterations(), and AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::Solve().
std::vector<c_vector<double,DIM> > AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mDeformedPosition [protected] |
Deformed position: mDeformedPosition[i](j) = x_j for node i
Definition at line 193 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::rGetDeformedPosition().
std::vector<c_vector<double,DIM> > AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mSurfaceTractions [protected] |
The surface tractions (which should really be non-zero) for the boundary elements in mBoundaryElements.
Definition at line 199 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetSurfaceTractionBoundaryConditions().
c_vector<double,DIM>(* AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, 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 AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetFunctionalBodyForce().
c_vector<double,DIM>(* AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mpTractionBoundaryConditionFunction)(c_vector< double, DIM > &X, double t) [protected] |
An optionally provided (pointer to a) function, giving the surface traction as a function of undeformed position.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetFunctionalTractionBoundaryCondition().
bool AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mUsingBodyForceFunction [protected] |
Whether the functional version of the body force is being used or not
Definition at line 211 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetFunctionalBodyForce().
bool AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mUsingTractionBoundaryConditionFunction [protected] |
Whether the functional version of the surface traction is being used or not
Definition at line 214 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::SetFunctionalTractionBoundaryCondition().
double AbstractNonlinearElasticitySolver< COMPRESSIBILITY_TYPE, DIM >::mCurrentTime [protected] |
This solver is for static problems, however the body force or surface tractions could be a function of time. The user should call SetCurrentTime() if this is the case
Reimplemented in AbstractCardiacMechanicsSolver< DIM >.
Definition at line 220 of file AbstractNonlinearElasticitySolver.hpp.
Referenced by AbstractNonlinearElasticitySolver< COMPRESSIBLE, DIM >::SetCurrentTime().