Chaste Commit::ca8ccdedf819b6e02855bc0e8e6f50bdecbc5208
|
#include <IncompressibleNonlinearElasticitySolver.hpp>
Public Member Functions | |
IncompressibleNonlinearElasticitySolver (AbstractTetrahedralMesh< DIM, DIM > &rQuadMesh, SolidMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory) | |
virtual | ~IncompressibleNonlinearElasticitySolver () |
Public Member Functions inherited from AbstractNonlinearElasticitySolver< DIM > | |
void | ComputeResidual (Vec currentGuess, Vec residualVector) |
void | ComputeJacobian (Vec currentGuess, Mat *pJacobian, Mat *pPreconditioner) |
AbstractNonlinearElasticitySolver (AbstractTetrahedralMesh< DIM, DIM > &rQuadMesh, SolidMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory, CompressibilityType compressibilityType) | |
virtual | ~AbstractNonlinearElasticitySolver () |
void | Solve (double tol=-1.0) |
void | SetIncludeActiveTension (bool includeActiveTension=true) |
unsigned | GetNumNewtonIterations () |
void | SetWriteOutputEachNewtonIteration (bool writeOutputEachNewtonIteration=true) |
void | SetKspAbsoluteTolerance (double kspAbsoluteTolerance) |
void | SetTakeFullFirstNewtonStep (bool takeFullFirstStep=true) |
void | SetUsePetscDirectSolve (bool usePetscDirectSolve=true) |
void | SetCurrentTime (double time) |
void | CreateCmguiOutput () |
void | WriteCurrentStrains (StrainType strainType, std::string fileName, int counterToAppend=-1) |
void | SetComputeAverageStressPerElementDuringSolve (bool setComputeAverageStressPerElement=true) |
void | WriteCurrentAverageElementStresses (std::string fileName, int counterToAppend=-1) |
std::vector< c_vector< double, DIM > > & | rGetSpatialSolution () |
std::vector< c_vector< double, DIM > > & | rGetDeformedPosition () |
c_matrix< double, DIM, DIM > | GetAverageStressPerElement (unsigned elementIndex) |
Public Member Functions inherited from AbstractContinuumMechanicsSolver< DIM > | |
AbstractContinuumMechanicsSolver (AbstractTetrahedralMesh< DIM, DIM > &rQuadMesh, ContinuumMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory, CompressibilityType compressibilityType) | |
virtual | ~AbstractContinuumMechanicsSolver () |
void | WriteCurrentSpatialSolution (std::string fileName, std::string fileExtension, int counterToAppend=-1) |
void | WriteCurrentPressureSolution (int counterToAppend=-1) |
void | SetWriteOutput (bool writeOutput=true) |
void | CreateVtkOutput (std::string spatialSolutionName="Spatial solution") |
std::vector< double > & | rGetCurrentSolution () |
std::vector< double > & | rGetPressures () |
Protected Member Functions | |
virtual void | AssembleOnElement (Element< DIM, DIM > &rElement, c_matrix< double, STENCIL_SIZE, STENCIL_SIZE > &rAElem, c_matrix< double, STENCIL_SIZE, STENCIL_SIZE > &rAElemPrecond, c_vector< double, STENCIL_SIZE > &rBElem, bool assembleResidual, bool assembleJacobian) |
void | FormInitialGuess () |
void | AssembleSystem (bool assembleResidual, bool assembleJacobian) |
Protected Member Functions inherited from AbstractNonlinearElasticitySolver< DIM > | |
void | AddStressToAverageStressPerElement (c_matrix< double, DIM, DIM > &rT, unsigned elementIndex) |
virtual void | SetKspSolverAndPcType (KSP solver) |
virtual void | FinishAssembleSystem (bool assembleResidual, bool assembleLinearSystem) |
void | GetElementCentroidStrain (StrainType strainType, Element< DIM, DIM > &rElement, c_matrix< double, DIM, DIM > &rDeformationGradient) |
virtual void | AddActiveStressAndStressDerivative (c_matrix< double, DIM, DIM > &rC, unsigned elementIndex, unsigned currentQuadPointGlobalIndex, c_matrix< double, DIM, DIM > &rT, FourthOrderTensor< DIM, DIM, DIM, DIM > &rDTdE, bool addToDTdE) |
virtual void | SetupChangeOfBasisMatrix (unsigned elementIndex, unsigned currentQuadPointGlobalIndex) |
void | AssembleOnBoundaryElement (BoundaryElement< DIM-1, DIM > &rBoundaryElement, c_matrix< double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE > &rAelem, c_vector< double, BOUNDARY_STENCIL_SIZE > &rBelem, bool assembleResidual, bool assembleJacobian, unsigned boundaryConditionIndex) |
bool | ShouldAssembleMatrixTermForPressureOnDeformedBc () |
void | AssembleOnBoundaryElementForPressureOnDeformedBc (BoundaryElement< DIM-1, DIM > &rBoundaryElement, c_matrix< double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE > &rAelem, c_vector< double, BOUNDARY_STENCIL_SIZE > &rBelem, bool assembleResidual, bool assembleJacobian, unsigned boundaryConditionIndex) |
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) |
void | SolveNonSnes (double tol=-1.0) |
Protected Member Functions inherited from AbstractContinuumMechanicsSolver< DIM > | |
void | AllocateMatrixMemory () |
void | ApplyDirichletBoundaryConditions (ApplyDirichletBcsType type, bool symmetricProblem) |
void | AddIdentityBlockForDummyPressureVariables (ApplyDirichletBcsType type) |
void | RemovePressureDummyValuesThroughLinearInterpolation () |
Static Protected Attributes | |
static const size_t | NUM_NODES_PER_ELEMENT = AbstractNonlinearElasticitySolver<DIM>::NUM_NODES_PER_ELEMENT |
static const size_t | NUM_VERTICES_PER_ELEMENT = AbstractNonlinearElasticitySolver<DIM>::NUM_VERTICES_PER_ELEMENT |
static const size_t | NUM_NODES_PER_BOUNDARY_ELEMENT = AbstractNonlinearElasticitySolver<DIM>::NUM_NODES_PER_BOUNDARY_ELEMENT |
static const size_t | STENCIL_SIZE = DIM*NUM_NODES_PER_ELEMENT + NUM_VERTICES_PER_ELEMENT |
static const size_t | BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT |
Static Protected Attributes inherited from AbstractNonlinearElasticitySolver< DIM > | |
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 const size_t | BOUNDARY_STENCIL_SIZE = DIM*NUM_NODES_PER_BOUNDARY_ELEMENT |
static double | MAX_NEWTON_ABS_TOL = 1e-7 |
static double | MIN_NEWTON_ABS_TOL = 1e-10 |
static double | NEWTON_REL_TOL = 1e-4 |
Friends | |
class | TestIncompressibleNonlinearElasticitySolver |
class | TestCompressibleNonlinearElasticitySolver |
class | TestNonlinearElasticityAdjointSolver |
class | AdaptiveNonlinearElasticityProblem |
Finite elasticity solver. Solves static *incompressible* nonlinear elasticity problems with arbitrary (incompressible) material laws and a body force.
Uses quadratic-linear bases (for displacement and pressure), and is therefore outside other assembler or solver hierarchy.
Definition at line 60 of file IncompressibleNonlinearElasticitySolver.hpp.
IncompressibleNonlinearElasticitySolver< DIM >::IncompressibleNonlinearElasticitySolver | ( | AbstractTetrahedralMesh< DIM, DIM > & | rQuadMesh, |
SolidMechanicsProblemDefinition< DIM > & | rProblemDefinition, | ||
std::string | outputDirectory | ||
) |
Constructor.
rQuadMesh | The quadratic mesh to solve on |
rProblemDefinition | an object defining in particular the body force and boundary conditions |
outputDirectory | The output directory |
Definition at line 615 of file IncompressibleNonlinearElasticitySolver.cpp.
References EXCEPTION, IncompressibleNonlinearElasticitySolver< DIM >::FormInitialGuess(), and SolidMechanicsProblemDefinition< DIM >::GetCompressibilityType().
|
inlinevirtual |
Destructor.
Definition at line 155 of file IncompressibleNonlinearElasticitySolver.hpp.
|
protectedvirtual |
Assemble residual or Jacobian on an element, using the current solution stored in mCurrrentSolution. The ordering assumed is (in 2d) rBElem = [u0 v0 u1 v1 .. u5 v5 p0 p1 p2]. (This ordering in used at the element level, but not in the global matrix/vector).
rElement | The element to assemble on. |
rAElem | The element's contribution to the LHS matrix is returned in this n by n matrix, where n is the no. of nodes in this element. There is no need to zero this matrix before calling. |
rAElemPrecond | The element's contribution to the matrix passed to PetSC in creating a preconditioner |
rBElem | The element's contribution to the RHS vector is returned in this vector of length n, the no. of nodes in this element. There is no need to zero this vector before calling. |
assembleResidual | A bool stating whether to assemble the residual vector. |
assembleJacobian | A bool stating whether to assemble the Jacobian matrix. |
Definition at line 210 of file IncompressibleNonlinearElasticitySolver.cpp.
References LinearBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), QuadraticBasisFunction< ELEMENT_DIM >::ComputeBasisFunctions(), AbstractMaterialLaw< DIM >::ComputeStressAndStressDerivative(), QuadraticBasisFunction< ELEMENT_DIM >::ComputeTransformedBasisFunctionDerivatives(), Determinant(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), Inverse(), NEVER_REACHED, and AbstractMaterialLaw< DIM >::SetChangeOfBasisMatrix().
|
protectedvirtual |
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).
assembleResidual | A bool stating whether to assemble the residual vector. |
assembleJacobian | A bool stating whether to assemble the Jacobian matrix. |
Implements AbstractNonlinearElasticitySolver< DIM >.
Definition at line 52 of file IncompressibleNonlinearElasticitySolver.cpp.
References PetscVecTools::Finalise(), PetscTools::GetMyRank(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetNodeGlobalIndex(), AbstractElement< ELEMENT_DIM, SPACE_DIM >::GetOwnership(), CommandLineArguments::Instance(), CommandLineArguments::OptionExists(), PetscMatTools::SwitchWriteMode(), PetscMatTools::Zero(), and PetscVecTools::Zero().
|
protected |
Set up the current guess to be the solution given no displacement.
The initial guess is zero for spatial variables, 0 for dummy pressure variables, and the appropriate choice of pressure such that the initial stress is zero (depends on material law)
In a homogeneous problem, all (the non-dummy) pressures are the same. In a heterogeneous problem, p for a given vertex is the zero-strain-pressure for ONE of the elements containing that vertex (which element containing the vertex is reached LAST). In this case the initial guess will be close but not exactly the solution given zero body force.
Definition at line 589 of file IncompressibleNonlinearElasticitySolver.cpp.
Referenced by IncompressibleNonlinearElasticitySolver< DIM >::IncompressibleNonlinearElasticitySolver().
|
friend |
Definition at line 65 of file IncompressibleNonlinearElasticitySolver.hpp.
|
friend |
Definition at line 63 of file IncompressibleNonlinearElasticitySolver.hpp.
|
friend |
Definition at line 62 of file IncompressibleNonlinearElasticitySolver.hpp.
|
friend |
Definition at line 64 of file IncompressibleNonlinearElasticitySolver.hpp.
|
staticprotected |
Boundary stencil size.
Definition at line 85 of file IncompressibleNonlinearElasticitySolver.hpp.
|
staticprotected |
Number of nodes per boundary element.
Definition at line 76 of file IncompressibleNonlinearElasticitySolver.hpp.
|
staticprotected |
Number of nodes per element.
Definition at line 70 of file IncompressibleNonlinearElasticitySolver.hpp.
|
staticprotected |
Number of vertices per element.
Definition at line 73 of file IncompressibleNonlinearElasticitySolver.hpp.
|
staticprotected |
Stencil size - number of unknowns per element (DIM*NUM_NODES_PER_ELEMENT displacement unknowns, NUM_VERTICES_PER_ELEMENT pressure unknowns.
Definition at line 82 of file IncompressibleNonlinearElasticitySolver.hpp.