Chaste Commit::ca8ccdedf819b6e02855bc0e8e6f50bdecbc5208
|
#include <AbstractContinuumMechanicsSolver.hpp>
Public Member Functions | |
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 () |
virtual std::vector< c_vector< double, DIM > > & | rGetSpatialSolution ()=0 |
std::vector< double > & | rGetPressures () |
Protected Member Functions | |
void | AllocateMatrixMemory () |
void | ApplyDirichletBoundaryConditions (ApplyDirichletBcsType type, bool symmetricProblem) |
void | AddIdentityBlockForDummyPressureVariables (ApplyDirichletBcsType type) |
void | RemovePressureDummyValuesThroughLinearInterpolation () |
Protected Attributes | |
AbstractTetrahedralMesh< DIM, DIM > & | mrQuadMesh |
ContinuumMechanicsProblemDefinition< DIM > & | mrProblemDefinition |
bool | mWriteOutput |
std::string | mOutputDirectory |
OutputFileHandler * | mpOutputFileHandler |
std::vector< c_vector< double, DIM > > | mSpatialSolution |
std::vector< double > | mPressureSolution |
std::vector< double > | mCurrentSolution |
GaussianQuadratureRule< DIM > * | mpQuadratureRule |
GaussianQuadratureRule< DIM-1 > * | mpBoundaryQuadratureRule |
CompressibilityType | mCompressibilityType |
unsigned | mProblemDimension |
unsigned | mNumDofs |
bool | mVerbose |
Vec | mResidualVector |
Vec | mLinearSystemRhsVector |
Mat | mSystemLhsMatrix |
Vec | mDirichletBoundaryConditionsVector |
Mat | mPreconditionMatrix |
Friends | |
class | StressRecoveror< DIM > |
class | VtkNonlinearElasticitySolutionWriter< DIM > |
General base class for continuum mechanics solvers. Deals with memory allocation, writing to file, applying boundary conditions, and adding an identity block (see below)
Note: the PROBLEM DIMENSION depends on DIM and whether a compressible or incompressible problem is being solved. (i) Compressible: mProblemDimension is set to DIM, and the ordering of the unknowns is, for 2D say: [U1 V1 U2 V2 .. Un Vn] (ii) Incompressible: Here we have spatial unknowns (at all nodes of the QuadraticMesh) and pressure unknowns (only at the vertices, as linear bases are). We choose mProblemDim = DIM+1, use (for parallelisation reasons) the ordering [U1 V1 P1 , .. , Un Vn, Pn], introducing dummy variables for pressure unknowns at internal hences, with the equation Pi=0 at these nodes, hence the need for an identity block to be added the corresponding parts of the matrix.
Definition at line 87 of file AbstractContinuumMechanicsSolver.hpp.
AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver | ( | AbstractTetrahedralMesh< DIM, DIM > & | rQuadMesh, |
ContinuumMechanicsProblemDefinition< DIM > & | rProblemDefinition, | ||
std::string | outputDirectory, | ||
CompressibilityType | compressibilityType | ||
) |
Constructor
rQuadMesh | the mesh |
rProblemDefinition | problem definition object |
outputDirectory | output directory name |
compressibilityType | 'INCOMPRESSIBLE' or 'COMPRESSIBLE' |
Definition at line 376 of file AbstractContinuumMechanicsSolver.hpp.
References AbstractContinuumMechanicsSolver< DIM >::AllocateMatrixMemory(), EXCEPTION, AbstractMesh< ELEMENT_DIM, SPACE_DIM >::GetNumNodes(), CommandLineArguments::Instance(), AbstractContinuumMechanicsSolver< DIM >::mCompressibilityType, AbstractContinuumMechanicsSolver< DIM >::mCurrentSolution, AbstractContinuumMechanicsSolver< DIM >::mNumDofs, AbstractContinuumMechanicsSolver< DIM >::mOutputDirectory, AbstractContinuumMechanicsSolver< DIM >::mpBoundaryQuadratureRule, AbstractContinuumMechanicsSolver< DIM >::mpOutputFileHandler, AbstractContinuumMechanicsSolver< DIM >::mpQuadratureRule, AbstractContinuumMechanicsSolver< DIM >::mProblemDimension, AbstractContinuumMechanicsSolver< DIM >::mrProblemDefinition, AbstractContinuumMechanicsSolver< DIM >::mrQuadMesh, AbstractContinuumMechanicsSolver< DIM >::mVerbose, AbstractContinuumMechanicsSolver< DIM >::mWriteOutput, and CommandLineArguments::OptionExists().
|
virtual |
Destructor
Definition at line 430 of file AbstractContinuumMechanicsSolver.hpp.
References PetscTools::Destroy().
|
protected |
For incompressible problems, we use the following ordering: [U0 V0 W0 P0 U1 V1 W1 P1 .. Un Vn Wn Pn] where (U V W) is the displacement, and P is the pressure. Therefore P is defined at every node; however for the solve, linear basis functions are used for the pressure, so pressure is solved for at each vertex, not at internal nodes (which are for quadratic basis functions).
The pressure variables for non-vertex nodes are therefore dummy variables. This method enforces the condition P_i=0, where i corresponds to a non-vertex node.
The first input parameter should be one of the following LINEAR_PROBLEM – indicating the overall problem is linear, and in which case the matrix will be altered (1 on diagonal, zeros assumed on rest of row) and the rhs vector will be set to 0.
NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY – indicating the overall problem is nonlinear and here only the residual vector will be altered (sets the residual value of the appropriate rows to P_i-0.0.
NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING – indicating the overall problem is nonlinear, and here, the residual vector will be altered, as will the matrix and RHS vector (see documentation on mResidualVector for why there is a separate residual vector and RHS vector).
Note the row (and columns) for the dummy variables are not explicitly zeroed, as they are assumed to be zero already - in a finite element assembly nothing will have been written to the rows or columns corresponding to dummy variables. Hence this method always maintains symmetry.
type | see above |
Definition at line 798 of file AbstractContinuumMechanicsSolver.hpp.
References PetscMatTools::SetElement(), and PetscVecTools::SetElement().
|
protected |
Allocates memory for the matrices and vectors
We want to allocate different numbers of non-zeros per row, which means PetscTools::SetupMat isn't that useful. We could call
but we would get warnings due to the lack allocation
Definition at line 830 of file AbstractContinuumMechanicsSolver.hpp.
References PetscTools::Destroy(), PetscTools::IsSequential(), and PetscTools::SetupMat().
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
Apply the Dirichlet boundary conditions to the linear system.
The first input parameter should be one of the following LINEAR_PROBLEM – indicating the overall problem is linear, and in which case the BCs will be applied to both the matrix and vector of the linear system
NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY – indicating the overall problem is nonlinear and here only the residual vector will be altered (apply the Dirichlet boundary conditions to the residual vector involves setting appropriate components to the difference between the current value and the correct value).
NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING – indicating the overall problem is nonlinear, and here, the residual vector will be altered, as will the matrix and RHS vector (see documentation on mResidualVector for why there is a separate residual vector and RHS vector).
The second parameter should be true if the overall problem is symmetric, in which case boundary conditions will be applied keeping the matrix symmetric (if the matrix is being altered). See in-code comments for how this is done.
type | see above |
symmetricProblem | see above |
Definition at line 672 of file AbstractContinuumMechanicsSolver.hpp.
References PetscVecTools::AddScaledVector(), PetscTools::Destroy(), PetscMatTools::Finalise(), PetscVecTools::Finalise(), PetscMatTools::GetMatrixRowDistributed(), PetscVecTools::SetElement(), PetscVecTools::Zero(), PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal(), and PetscMatTools::ZeroRowsWithValueOnDiagonal().
void AbstractContinuumMechanicsSolver< DIM >::CreateVtkOutput | ( | std::string | spatialSolutionName = "Spatial solution" | ) |
Convert the output to vtk format (placed in a folder called vtk in the output directory).
spatialSolutionName | is used to identify the spatial solution as a velocity, displacement... |
** TO BE DEPRECATED - see #2321 **
Definition at line 547 of file AbstractContinuumMechanicsSolver.hpp.
References VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::AddCellData(), VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::AddPointData(), EXCEPTION, and VtkMeshWriter< ELEMENT_DIM, SPACE_DIM >::WriteFilesUsingMesh().
|
protected |
For incompressible problems, we use the following ordering: [U0 V0 W0 P0 U1 V1 W1 P1 .. Un Vn Wn Pn] where (U V W) is the displacement, and P is the pressure. Therefore P is defined at every node; however for the solve, linear basis functions are used for the pressure, so pressure is solved for at each vertex, not at internal nodes (which are for quadratic basis functions).
The above method, AddIdentityBlockForDummyPressureVariables(), enforces the condition P_i=0, where i corresponds to a non-vertex node.
AFTER the solve, this method can be used to remove the dummy values, but looping over all edges, and linearly interpolating the pressure at the two vertices onto the internal node.
This method assumes each internal node is midway between the two vertices.
Definition at line 593 of file AbstractContinuumMechanicsSolver.hpp.
|
inline |
Definition at line 353 of file AbstractContinuumMechanicsSolver.hpp.
References AbstractContinuumMechanicsSolver< DIM >::mCurrentSolution.
std::vector< double > & AbstractContinuumMechanicsSolver< DIM >::rGetPressures | ( | ) |
Only valid if mCompressibilityType==INCOMPRESSIBLE
Definition at line 578 of file AbstractContinuumMechanicsSolver.hpp.
|
pure virtual |
Implemented in AbstractNonlinearElasticitySolver< DIM >, and StokesFlowSolver< DIM >.
void AbstractContinuumMechanicsSolver< DIM >::SetWriteOutput | ( | bool | writeOutput = true | ) |
Set whether to write any output.
writeOutput | (defaults to true) |
Definition at line 536 of file AbstractContinuumMechanicsSolver.hpp.
References EXCEPTION.
void AbstractContinuumMechanicsSolver< DIM >::WriteCurrentPressureSolution | ( | int | counterToAppend = -1 | ) |
Write the pressure solution. Only valid if mCompressibilityType==INCOMPRESSIBLE. Writes the pressure for ALL nodes on the mesh, including internal nodes (as these are not assumed to be have indices greater than vertex nodes). As linear basis functions are used for pressure, the pressure solution is only computed at the vertices, and pressure dummy variables are used at internal nodes, hence for each internal node, 0 will be written to file.
counterToAppend | append a counter to the file name |
WriteCurrentPressureSolution() –> file called "pressure.txt" WriteCurrentPressureSolution(3) –> file called "pressure_3.txt"
Definition at line 500 of file AbstractContinuumMechanicsSolver.hpp.
References PetscTools::AmMaster(), and PetscTools::Barrier().
void AbstractContinuumMechanicsSolver< DIM >::WriteCurrentSpatialSolution | ( | std::string | fileName, |
std::string | fileExtension, | ||
int | counterToAppend = -1 |
||
) |
Write the spatial solution (deformed position if solids, flow if fluids) at the nodes
fileName | (stem) |
fileExtension | to append at end. |
counterToAppend | append a counter to the file name (defaults to nothing appended). |
For example: WriteCurrentSpatialSolution("solution","nodes") –> file called "solution.nodes" WriteCurrentSpatialSolution("solution","nodes",3) –> file called "solution_3.nodes"
Definition at line 458 of file AbstractContinuumMechanicsSolver.hpp.
References PetscTools::AmMaster(), and PetscTools::Barrier().
|
friend |
Definition at line 830 of file AbstractContinuumMechanicsSolver.hpp.
|
friend |
Definition at line 830 of file AbstractContinuumMechanicsSolver.hpp.
|
protected |
This is equal to either COMPRESSIBLE or INCOMPRESSIBLE (see enumeration class) and is only used in computing mNumDofs and allocating matrix memory.
Definition at line 138 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
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 126 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver(), and AbstractContinuumMechanicsSolver< DIM >::rGetCurrentSolution().
|
protected |
Helper vector (see ApplyDirichletBoundaryConditions code).
Definition at line 197 of file AbstractContinuumMechanicsSolver.hpp.
|
protected |
The RHS side in the linear system that is solved each Newton iteration.
Definition at line 187 of file AbstractContinuumMechanicsSolver.hpp.
|
protected |
Number of degrees of freedom (equal to mProblemDim*num_nodes)
Definition at line 153 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
Where to write output, relative to CHASTE_TEST_OUTPUT.
Definition at line 106 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
Boundary Gaussian quadrature rule.
Definition at line 132 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
Output file handler.
Definition at line 109 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
Gaussian quadrature rule.
Definition at line 129 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
Precondition matrix for the linear system.
Definition at line 202 of file AbstractContinuumMechanicsSolver.hpp.
|
protected |
Pressures solution at each vertex of the mesh. Only valid if mCompressibilityType==INCOMPRESSIBLE.
Definition at line 117 of file AbstractContinuumMechanicsSolver.hpp.
|
protected |
The problem dimension - the number of unknowns at each node. For compressible problems, where only the deformation/flow is computed at each node, this is equal to DIM For incompressible problems, where pressure is computed as well, this is equal to DIM+1 (there is a pressure variable defined even for internal nodes at which pressure is not computed, this is a dummy pressure variable – done like this for parallelisation reasons. See above).
Definition at line 148 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
Residual vector nonlinear problems.
Since the residual in nonlinear problems is usually also the RHS vector in the linear system, it may seem unncessary to also have the member variable mLinearSystemRhsVector.
However: Newton's method is Ju = f, where J is the Jacobian, u the (negative of the) update and f the residual, but when applying Dirichlet boundary conditions in the compressible case, we alter the rows of the matrix and 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
mLinearSystemRhsVector below represents f*.
Definition at line 182 of file AbstractContinuumMechanicsSolver.hpp.
|
protected |
Problem definition class - contains info on boundary conditions, etc
Definition at line 100 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
The mesh to be solved on. Requires 6 nodes per triangle (or 10 per tetrahedron) as quadratic bases are used.
Definition at line 97 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver(), and StokesFlowSolver< DIM >::StokesFlowSolver().
|
protected |
Spatial solution: For solids problems, mSpatialSolution[i](j) = x_j (new position) for node i. For fluids problems, mSpatialSolution[i](j) = u_j (flow) for node i.
Definition at line 114 of file AbstractContinuumMechanicsSolver.hpp.
|
protected |
Jacobian matrix of the nonlinear system, LHS matrix for the linear system.
Definition at line 192 of file AbstractContinuumMechanicsSolver.hpp.
|
protected |
If SetVerboseDuringSolve() is called on the problem definition class, or the command line argument "-mech_verbose" or "-mech_very_verbose" is given, than this bool will be set to true and lots of details about each nonlinear solve (including timing breakdowns) will be printed out
Definition at line 160 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().
|
protected |
Whether to write any output.
Definition at line 103 of file AbstractContinuumMechanicsSolver.hpp.
Referenced by AbstractContinuumMechanicsSolver< DIM >::AbstractContinuumMechanicsSolver().