Chaste
Release::3.4
|
#include <LinearSystem.hpp>
Public Member Functions | |
LinearSystem (PetscInt lhsVectorSize, unsigned rowPreallocation=UINT_MAX) | |
LinearSystem (Vec templateVector, unsigned rowPreallocation, bool newAllocationError=true) | |
LinearSystem (Vec residualVector, Mat jacobianMatrix) | |
LinearSystem (PetscInt lhsVectorSize, Mat lhsMatrix, Vec rhsVector) | |
~LinearSystem () | |
void | SetMatrixElement (PetscInt row, PetscInt col, double value) |
void | AddToMatrixElement (PetscInt row, PetscInt col, double value) |
void | AssembleFinalLinearSystem () |
void | AssembleIntermediateLinearSystem () |
void | FinaliseLhsMatrix () |
void | FinalisePrecondMatrix () |
void | SwitchWriteModeLhsMatrix () |
void | FinaliseRhsVector () |
void | SetMatrixIsSymmetric (bool isSymmetric=true) |
bool | IsMatrixSymmetric () |
void | SetMatrixIsConstant (bool matrixIsConstant) |
void | SetRelativeTolerance (double relativeTolerance) |
void | SetAbsoluteTolerance (double absoluteTolerance) |
void | SetKspType (const char *kspType) |
void | SetPcType (const char *pcType, boost::shared_ptr< std::vector< PetscInt > > pBathNodes=boost::shared_ptr< std::vector< PetscInt > >()) |
void | DisplayMatrix () |
void | DisplayRhs () |
void | SetMatrixRow (PetscInt row, double value) |
Vec | GetMatrixRowDistributed (unsigned rowIndex) |
void | ZeroMatrixRowsWithValueOnDiagonal (std::vector< unsigned > &rRows, double diagonalValue) |
void | ZeroMatrixRowsAndColumnsWithValueOnDiagonal (std::vector< unsigned > &rRowColIndices, double diagonalValue) |
void | ZeroMatrixColumn (PetscInt col) |
void | ZeroLhsMatrix () |
void | ZeroRhsVector () |
void | ZeroLinearSystem () |
Vec | Solve (Vec lhsGuess=NULL) |
void | SetRhsVectorElement (PetscInt row, double value) |
void | AddToRhsVectorElement (PetscInt row, double value) |
unsigned | GetSize () const |
void | SetNullBasis (Vec nullbasis[], unsigned numberOfBases) |
void | RemoveNullSpace () |
Vec & | rGetRhsVector () |
Vec | GetRhsVector () const |
Mat & | rGetLhsMatrix () |
Mat | GetLhsMatrix () const |
Mat & | rGetPrecondMatrix () |
Vec & | rGetDirichletBoundaryConditionsVector () |
void | GetOwnershipRange (PetscInt &lo, PetscInt &hi) |
double | GetMatrixElement (PetscInt row, PetscInt col) |
double | GetRhsVectorElement (PetscInt row) |
unsigned | GetNumIterations () const |
template<size_t MATRIX_SIZE> | |
void | AddLhsMultipleValues (unsigned *matrixRowAndColIndices, c_matrix< double, MATRIX_SIZE, MATRIX_SIZE > &rSmallMatrix) |
template<size_t VECTOR_SIZE> | |
void | AddRhsMultipleValues (unsigned *vectorIndices, c_vector< double, VECTOR_SIZE > &smallVector) |
void | SetPrecondMatrixIsDifferentFromLhs (bool precondIsDifferent=true) |
void | SetUseFixedNumberIterations (bool useFixedNumberIterations=true, unsigned evaluateNumItsEveryNSolves=UINT_MAX) |
void | ResetKspSolver () |
Private Member Functions | |
template<class Archive > | |
void | serialize (Archive &archive, const unsigned int version) |
Private Attributes | |
Mat | mLhsMatrix |
Mat | mPrecondMatrix |
Vec | mRhsVector |
PetscInt | mSize |
PetscInt | mOwnershipRangeLo |
PetscInt | mOwnershipRangeHi |
MatNullSpace | mMatNullSpace |
bool | mDestroyMatAndVec |
KSP | mKspSolver |
bool | mKspIsSetup |
double | mNonZerosUsed |
bool | mMatrixIsConstant |
double | mTolerance |
bool | mUseAbsoluteTolerance |
std::string | mKspType |
std::string | mPcType |
Vec | mDirichletBoundaryConditionsVector |
PCBlockDiagonal * | mpBlockDiagonalPC |
PCLDUFactorisation * | mpLDUFactorisationPC |
PCTwoLevelsBlockDiagonal * | mpTwoLevelsBlockDiagonalPC |
boost::shared_ptr< std::vector < PetscInt > > | mpBathNodes |
bool | mPrecondMatrixIsNotLhs |
unsigned | mRowPreallocation |
bool | mUseFixedNumberIterations |
unsigned | mEvaluateNumItsEveryNSolves |
void * | mpConvergenceTestContext |
unsigned | mNumSolves |
PetscReal | mEigMin |
PetscReal | mEigMax |
bool | mForceSpectrumReevaluation |
Friends | |
class | TestLinearSystem |
class | TestPCBlockDiagonal |
class | TestPCTwoLevelsBlockDiagonal |
class | TestPCLDUFactorisation |
class | TestChebyshevIteration |
class | boost::serialization::access |
Linear System class. Stores and solves a linear equation of the form Ax=b, where A is a square matrix and x and b are column vectors. The class uses PETSc.
Definition at line 65 of file LinearSystem.hpp.
Constructor.
lhsVectorSize | the size of the LHS vector |
rowPreallocation | the max number of nonzero entries expected on a row
|
Definition at line 54 of file LinearSystem.cpp.
References PetscTools::CreateVec(), EXCEPTION, mKspType, mLhsMatrix, mNumSolves, mOwnershipRangeHi, mOwnershipRangeLo, mPcType, mRhsVector, mRowPreallocation, mSize, and PetscTools::SetupMat().
LinearSystem::LinearSystem | ( | Vec | templateVector, |
unsigned | rowPreallocation, | ||
bool | newAllocationError = true |
||
) |
Alternative constructor.
Create a linear system, where the size is based on the size of a given PETSc vec.
The LHS & RHS vectors will be created by duplicating this vector's settings. This should avoid problems with using VecScatter on bidomain simulation results.
templateVector | a PETSc vec |
rowPreallocation | the max number of nonzero entries expected on a row |
newAllocationError | tells PETSc whether to set the MAT_NEW_NONZERO_ALLOCATION_ERR. ** currently only used in PETSc 3.3 and later ** in PETSc 3.2 and earlier MAT_NEW_NONZERO_ALLOCATION_ERR defaults to false in PETSc 3.3 MAT_NEW_NONZERO_ALLOCATION_ERR defaults to true |
Definition at line 146 of file LinearSystem.cpp.
References mKspType, mLhsMatrix, mNumSolves, mOwnershipRangeHi, mOwnershipRangeLo, mPcType, mRhsVector, mRowPreallocation, mSize, and PetscTools::SetupMat().
Alternative constructor.
Create a linear system which wraps the provided PETSc objects so we can access them using our API. Either of the objects may be NULL, but at least one of them must not be.
Useful for storing residuals and jacobians when solving nonlinear PDEs.
residualVector | the residual vector |
jacobianMatrix | the Jacobian matrix |
Definition at line 187 of file LinearSystem.cpp.
References mKspType, mLhsMatrix, mNumSolves, mOwnershipRangeHi, mOwnershipRangeLo, mPcType, mRhsVector, mRowPreallocation, and mSize.
Alternative constructor for archiving.
lhsVectorSize | the size of the LHS vector |
lhsMatrix | the RHS matrix |
rhsVector | the RHS vector |
Definition at line 109 of file LinearSystem.cpp.
References mLhsMatrix, mNumSolves, mOwnershipRangeHi, mOwnershipRangeLo, and mRhsVector.
LinearSystem::~LinearSystem | ( | ) |
Destructor.
Definition at line 251 of file LinearSystem.cpp.
References PetscTools::AmMaster(), PetscTools::Destroy(), mDestroyMatAndVec, mDirichletBoundaryConditionsVector, mKspIsSetup, mKspSolver, mLhsMatrix, mMatNullSpace, mNumSolves, mpBlockDiagonalPC, mpConvergenceTestContext, mpLDUFactorisationPC, mPrecondMatrix, mPrecondMatrixIsNotLhs, mpTwoLevelsBlockDiagonalPC, mRhsVector, and PETSC_DESTROY_PARAM.
|
inline |
Add multiple values to the matrix of linear system.
matrixRowAndColIndices | mapping from index of the ublas matrix (see param below) to index of the PETSc matrix of this linear system |
rSmallMatrix | Ublas matrix containing the values to be added |
N.B. Values which are not local (ie the row is not owned) will be skipped.
Definition at line 540 of file LinearSystem.hpp.
References PetscMatTools::AddMultipleValues(), and mLhsMatrix.
|
inline |
Add multiple values to the RHS vector.
vectorIndices | mapping from index of the ublas vector (see param below) to index of the vector of this linear system |
smallVector | Ublas vector containing the values to be added |
N.B. Values which are not local (ie the row is not owned) will be skipped.
Definition at line 555 of file LinearSystem.hpp.
References PetscVecTools::AddMultipleValues(), and mRhsVector.
Add the specified value to an entry of the matrix.
row | the row index |
col | the column index |
value | the value for this entry |
Definition at line 315 of file LinearSystem.cpp.
References PetscMatTools::AddToElement(), and mLhsMatrix.
Add a value to an element of the right-hand side vector.
row | the row index |
value | the value to set this entry |
Definition at line 357 of file LinearSystem.cpp.
References PetscVecTools::AddToElement(), and mRhsVector.
void LinearSystem::AssembleFinalLinearSystem | ( | ) |
Call this before Solve().
This calls FinaliseLhsMatrix() and FinaliseRhsVector().
Definition at line 320 of file LinearSystem.cpp.
References FinaliseLhsMatrix(), and FinaliseRhsVector().
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().
void LinearSystem::AssembleIntermediateLinearSystem | ( | ) |
Should be called before AddToMatrixElement.
This calls SwitchWriteModeLhsMatrix() and FinaliseRhsVector().
Definition at line 326 of file LinearSystem.cpp.
References FinaliseRhsVector(), and SwitchWriteModeLhsMatrix().
void LinearSystem::DisplayMatrix | ( | ) |
Display the left-hand side matrix.
Definition at line 362 of file LinearSystem.cpp.
References PetscMatTools::Display(), and mLhsMatrix.
void LinearSystem::DisplayRhs | ( | ) |
Display the right-hand side vector.
Definition at line 367 of file LinearSystem.cpp.
References PetscVecTools::Display(), and mRhsVector.
void LinearSystem::FinaliseLhsMatrix | ( | ) |
Sets up the PETSc matrix left-hand-side mLhsMatrix
Definition at line 332 of file LinearSystem.cpp.
References PetscMatTools::Finalise(), and mLhsMatrix.
Referenced by AssembleFinalLinearSystem(), and AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::SetupGivenLinearSystem().
void LinearSystem::FinalisePrecondMatrix | ( | ) |
Sets up the PETSc matrix used for preconditioning.
Definition at line 342 of file LinearSystem.cpp.
References PetscMatTools::Finalise(), and mPrecondMatrix.
void LinearSystem::FinaliseRhsVector | ( | ) |
Sets up the PETSc vector right-hand-side mRhsVector
Definition at line 347 of file LinearSystem.cpp.
References PetscVecTools::Finalise(), and mRhsVector.
Referenced by AssembleFinalLinearSystem(), AssembleIntermediateLinearSystem(), and AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::SetupGivenLinearSystem().
Mat LinearSystem::GetLhsMatrix | ( | ) | const |
Definition at line 524 of file LinearSystem.cpp.
References mLhsMatrix.
row | the row index |
col | the column index |
Definition at line 489 of file LinearSystem.cpp.
References PetscMatTools::GetElement(), and mLhsMatrix.
Returns the i-th row of the LHS matrix as a distributed PETSc Vec
rowIndex | the row index |
Definition at line 377 of file LinearSystem.cpp.
References PetscMatTools::GetMatrixRowDistributed(), and mLhsMatrix.
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().
unsigned LinearSystem::GetNumIterations | ( | ) | const |
Definition at line 499 of file LinearSystem.cpp.
References mKspIsSetup, and mKspSolver.
lo | lowest index owned by this process |
hi | highest index owned by this process |
Definition at line 483 of file LinearSystem.cpp.
References mOwnershipRangeHi, and mOwnershipRangeLo.
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().
Vec LinearSystem::GetRhsVector | ( | ) | const |
Definition at line 514 of file LinearSystem.cpp.
References mRhsVector.
row | the row index |
Definition at line 494 of file LinearSystem.cpp.
References PetscVecTools::GetElement(), and mRhsVector.
unsigned LinearSystem::GetSize | ( | ) | const |
Definition at line 413 of file LinearSystem.cpp.
References mSize.
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().
bool LinearSystem::IsMatrixSymmetric | ( | ) |
Get whether PETSc considers the matrix in this linear system as symmetric or not.
Definition at line 568 of file LinearSystem.cpp.
References mLhsMatrix.
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().
void LinearSystem::RemoveNullSpace | ( | ) |
Remove the null space from the linear system.
Use for example if Dirichlet BC are applied to a singular system and, therefore, there's no null space anymore.
Definition at line 457 of file LinearSystem.cpp.
References mKspIsSetup, mKspSolver, mLhsMatrix, mMatNullSpace, PETSC_DESTROY_PARAM, and ResetKspSolver().
void LinearSystem::ResetKspSolver | ( | ) |
Method to regenerate all KSP objects, including the solver and the preconditioner (e.g. after changing the PDE time step when using time adaptivity).
Definition at line 1236 of file LinearSystem.cpp.
References mForceSpectrumReevaluation, mKspIsSetup, mKspSolver, and PETSC_DESTROY_PARAM.
Referenced by RemoveNullSpace().
Vec & LinearSystem::rGetDirichletBoundaryConditionsVector | ( | ) |
Should only be used by the BoundaryConditionsContainer.
Definition at line 539 of file LinearSystem.cpp.
References mDirichletBoundaryConditionsVector.
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().
Mat & LinearSystem::rGetLhsMatrix | ( | ) |
Definition at line 519 of file LinearSystem.cpp.
References mLhsMatrix.
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyPeriodicBcsToLinearProblem(), OdeLinearSystemSolver::rGetLhsMatrix(), AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::SetupGivenLinearSystem(), SimpleNewtonNonlinearSolver::Solve(), and OdeLinearSystemSolver::SolveOneTimeStep().
Mat & LinearSystem::rGetPrecondMatrix | ( | ) |
Definition at line 529 of file LinearSystem.cpp.
References EXCEPTION, mPrecondMatrix, and mPrecondMatrixIsNotLhs.
Vec & LinearSystem::rGetRhsVector | ( | ) |
Definition at line 509 of file LinearSystem.cpp.
References mRhsVector.
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem(), AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::SetupGivenLinearSystem(), SimpleNewtonNonlinearSolver::Solve(), and OdeLinearSystemSolver::SolveOneTimeStep().
|
inlineprivate |
Archive the member variables.
archive | |
version |
Definition at line 162 of file LinearSystem.hpp.
References mKspType, mMatrixIsConstant, mNonZerosUsed, mPcType, mTolerance, and mUseAbsoluteTolerance.
void LinearSystem::SetAbsoluteTolerance | ( | double | absoluteTolerance | ) |
Set the absolute tolerance.
absoluteTolerance | the absolute tolerance |
Definition at line 594 of file LinearSystem.cpp.
References mKspIsSetup, mKspSolver, mTolerance, and mUseAbsoluteTolerance.
void LinearSystem::SetKspType | ( | const char * | kspType | ) |
Set the KSP solver type (see PETSc KSPSetType() for valid arguments).
kspType | the KSP solver type |
Definition at line 604 of file LinearSystem.cpp.
References mKspIsSetup, mKspSolver, and mKspType.
Change one of the entires of the matrix to the specified value.
row | the row index |
col | the column index |
value | the value for this entry |
Definition at line 310 of file LinearSystem.cpp.
References mLhsMatrix, and PetscMatTools::SetElement().
void LinearSystem::SetMatrixIsConstant | ( | bool | matrixIsConstant | ) |
Set mMatrixIsConstant.
matrixIsConstant | whether the matrix is constant |
Definition at line 579 of file LinearSystem.cpp.
References mMatrixIsConstant.
void LinearSystem::SetMatrixIsSymmetric | ( | bool | isSymmetric = true | ) |
Force PETSc to treat the matrix in this linear system as symmetric from now on.
isSymmetric | whether the matrix is symmetric or not (defaults to true) |
Definition at line 544 of file LinearSystem.cpp.
References mLhsMatrix, and PetscMatTools::SetOption().
Set all entries in a given row of a matrix to a certain value. This must be called by the process who owns the row, (but other processors will treat it as a null-op)
row | the row index |
value | the value to set each entry in this row |
Definition at line 372 of file LinearSystem.cpp.
References mLhsMatrix, and PetscMatTools::SetRow().
Set the null basis of the linear system. In debug mode we test for orthonormality and throw EXCEPTIONs:
nullbasis | an array PETSc vectors containing orthogonal directions in the nullspace |
numberOfBases | the number of directions (size of nullbasis array) |
Definition at line 418 of file LinearSystem.cpp.
References EXCEPTION, and mMatNullSpace.
void LinearSystem::SetPcType | ( | const char * | pcType, |
boost::shared_ptr< std::vector< PetscInt > > | pBathNodes = boost::shared_ptr<std::vector<PetscInt> >() |
||
) |
Set the preconditioner type (see PETSc PCSetType() for valid arguments).
pcType | the preconditioner type |
pBathNodes | the list of nodes defining the bath |
Definition at line 614 of file LinearSystem.cpp.
References mKspIsSetup, mKspSolver, mpBathNodes, mpBlockDiagonalPC, mPcType, mpLDUFactorisationPC, mpTwoLevelsBlockDiagonalPC, and TERMINATE.
void LinearSystem::SetPrecondMatrixIsDifferentFromLhs | ( | bool | precondIsDifferent = true | ) |
Set method for mPrecondMatrixIsNotLhs.
precondIsDifferent | whether the matrix used for preconditioning is the same as the LHS. |
Definition at line 1204 of file LinearSystem.cpp.
References mOwnershipRangeHi, mOwnershipRangeLo, mPrecondMatrix, mPrecondMatrixIsNotLhs, mRowPreallocation, mSize, NEVER_REACHED, and PetscTools::SetupMat().
void LinearSystem::SetRelativeTolerance | ( | double | relativeTolerance | ) |
Set the relative tolerance.
relativeTolerance | the relative tolerance |
Definition at line 584 of file LinearSystem.cpp.
References mKspIsSetup, mKspSolver, mTolerance, and mUseAbsoluteTolerance.
Set an element of the right-hand side vector to a given value.
row | the row index |
value | the value to set this entry |
Definition at line 352 of file LinearSystem.cpp.
References mRhsVector, and PetscVecTools::SetElement().
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem(), and BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyPeriodicBcsToLinearProblem().
void LinearSystem::SetUseFixedNumberIterations | ( | bool | useFixedNumberIterations = true , |
unsigned | evaluateNumItsEveryNSolves = UINT_MAX |
||
) |
Set method for mUseFixedNumberIterations
useFixedNumberIterations | whether to use fixed number of iterations |
evaluateNumItsEveryNSolves | tells LinearSystem to perform a solve with convergence-based stop criteria every n solves to decide how many iterations perform for the next n-1 solves. Default is perfoming a single evaluation at the beginning of the simulation. |
Definition at line 1229 of file LinearSystem.cpp.
References mEvaluateNumItsEveryNSolves, and mUseFixedNumberIterations.
Solve the linear system.
lhsGuess | an optional initial guess for the solution (defaults to NULL) |
Definition at line 674 of file LinearSystem.cpp.
References PetscTools::AmMaster(), GenericEventHandler< 16, HeartEventHandler >::BeginEvent(), PetscTools::Destroy(), GenericEventHandler< 16, HeartEventHandler >::EndEvent(), EXCEPTION, mEigMax, mEigMin, mEvaluateNumItsEveryNSolves, mForceSpectrumReevaluation, mKspIsSetup, mKspSolver, mKspType, mLhsMatrix, mMatNullSpace, mMatrixIsConstant, mNonZerosUsed, mNumSolves, mpBathNodes, mpBlockDiagonalPC, mpConvergenceTestContext, mPcType, mpLDUFactorisationPC, mPrecondMatrix, mPrecondMatrixIsNotLhs, mpTwoLevelsBlockDiagonalPC, mRhsVector, mSize, mTolerance, mUseAbsoluteTolerance, mUseFixedNumberIterations, Timer::Print(), Timer::Reset(), TERMINATE, and PetscVecTools::Zero().
Referenced by SimpleNewtonNonlinearSolver::Solve(), and OdeLinearSystemSolver::SolveOneTimeStep().
void LinearSystem::SwitchWriteModeLhsMatrix | ( | ) |
Sets up the PETSc matrix left-hand-side mLhsMatrix
Definition at line 337 of file LinearSystem.cpp.
References mLhsMatrix, and PetscMatTools::SwitchWriteMode().
Referenced by AssembleIntermediateLinearSystem(), and AbstractAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, INTERPOLATION_LEVEL >::SetupGivenLinearSystem().
void LinearSystem::ZeroLhsMatrix | ( | ) |
Zero all entries of the left-hand side matrix.
Definition at line 402 of file LinearSystem.cpp.
References mLhsMatrix, and PetscMatTools::Zero().
Referenced by ZeroLinearSystem().
void LinearSystem::ZeroLinearSystem | ( | ) |
Zero all entries of the left-hand side matrix and right-hand side vector.
Definition at line 407 of file LinearSystem.cpp.
References ZeroLhsMatrix(), and ZeroRhsVector().
Referenced by SimpleNewtonNonlinearSolver::Solve().
void LinearSystem::ZeroMatrixColumn | ( | PetscInt | col | ) |
Zero a column of the left-hand side matrix.
Unfortunately there is no equivalent method in Petsc, so this has to be done carefully to ensure that the sparsity structure of the matrix is not broken. Only owned entries which are non-zero are zeroed.
col | the column index |
Definition at line 392 of file LinearSystem.cpp.
References mLhsMatrix, and PetscMatTools::ZeroColumn().
void LinearSystem::ZeroMatrixRowsAndColumnsWithValueOnDiagonal | ( | std::vector< unsigned > & | rRowColIndices, |
double | diagonalValue | ||
) |
Zero several rows and columns of the matrix, putting a given value on the diagonal.
rRowColIndices | A list of indices. All the rows with these indices, and all the columns with these indices, will be zeroed. |
diagonalValue | value to put in the diagonal entries (of the zeroed rows) |
Definition at line 387 of file LinearSystem.cpp.
References mLhsMatrix, and PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal().
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().
void LinearSystem::ZeroMatrixRowsWithValueOnDiagonal | ( | std::vector< unsigned > & | rRows, |
double | diagonalValue | ||
) |
Zero several rows of the matrix, putting a given value in the diagonal entries.
*Massively* less expensive than zeroing each matrix row individually
rRows | std::vector of rows to be zeroed |
diagonalValue | value to put in the diagonal entries (of the zeroed rows) |
Definition at line 382 of file LinearSystem.cpp.
References mLhsMatrix, and PetscMatTools::ZeroRowsWithValueOnDiagonal().
Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem(), and BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyPeriodicBcsToLinearProblem().
void LinearSystem::ZeroRhsVector | ( | ) |
Zero all entries of the right-hand side vector.
Definition at line 397 of file LinearSystem.cpp.
References mRhsVector, and PetscVecTools::Zero().
Referenced by ZeroLinearSystem().
|
friend |
Needed for serialization.
Definition at line 154 of file LinearSystem.hpp.
|
private |
Whether we need to destroy the PETSc matrix and vector in our destructor
Definition at line 90 of file LinearSystem.hpp.
Referenced by ~LinearSystem().
|
private |
Storage for efficient application of Dirichlet BCs, see AbstractBoundaryConditionsContainer
Definition at line 105 of file LinearSystem.hpp.
Referenced by rGetDirichletBoundaryConditionsVector(), and ~LinearSystem().
|
private |
Preconditioned operator largest eigenvalue
Definition at line 144 of file LinearSystem.hpp.
Referenced by Solve().
|
private |
Preconditioned operator smallest eigenvalue
Definition at line 141 of file LinearSystem.hpp.
Referenced by Solve().
|
private |
When using fixed number of iterations, a solve with residual-based stop criteria will be performed every mEvaluateNumItsEveryNSolves solves to decide how many iterations perform for the next mEvaluateNumItsEveryNSolves-1 solves
Definition at line 132 of file LinearSystem.hpp.
Referenced by SetUseFixedNumberIterations(), and Solve().
|
private |
Under certain circunstances you have to reevaluate the spectrum before the k*n-th, k=0,1,..., iteration
Definition at line 147 of file LinearSystem.hpp.
Referenced by ResetKspSolver(), and Solve().
|
private |
Used by Solve method to track whether KSP has been used.
Definition at line 93 of file LinearSystem.hpp.
Referenced by GetNumIterations(), RemoveNullSpace(), ResetKspSolver(), SetAbsoluteTolerance(), SetKspType(), SetPcType(), SetRelativeTolerance(), Solve(), and ~LinearSystem().
|
private |
The PETSc linear solver object
Definition at line 92 of file LinearSystem.hpp.
Referenced by GetNumIterations(), RemoveNullSpace(), ResetKspSolver(), SetAbsoluteTolerance(), SetKspType(), SetPcType(), SetRelativeTolerance(), Solve(), and ~LinearSystem().
|
private |
KSP solver type (see PETSc KSPSetType() )
Definition at line 102 of file LinearSystem.hpp.
Referenced by LinearSystem(), serialize(), SetKspType(), and Solve().
|
private |
The left-hand side matrix.
Definition at line 75 of file LinearSystem.hpp.
Referenced by AddLhsMultipleValues(), AddToMatrixElement(), DisplayMatrix(), FinaliseLhsMatrix(), GetLhsMatrix(), GetMatrixElement(), GetMatrixRowDistributed(), IsMatrixSymmetric(), LinearSystem(), RemoveNullSpace(), rGetLhsMatrix(), SetMatrixElement(), SetMatrixIsSymmetric(), SetMatrixRow(), Solve(), SwitchWriteModeLhsMatrix(), ZeroLhsMatrix(), ZeroMatrixColumn(), ZeroMatrixRowsAndColumnsWithValueOnDiagonal(), ZeroMatrixRowsWithValueOnDiagonal(), and ~LinearSystem().
|
private |
PETSc null matrix.
Definition at line 87 of file LinearSystem.hpp.
Referenced by RemoveNullSpace(), SetNullBasis(), Solve(), and ~LinearSystem().
|
private |
Whether the matrix is unchanged each time Solve() is called
Definition at line 95 of file LinearSystem.hpp.
Referenced by serialize(), SetMatrixIsConstant(), and Solve().
|
private |
Yes, it really is stored as a double.
Definition at line 94 of file LinearSystem.hpp.
Referenced by serialize(), and Solve().
|
private |
Number of solves performed since the current object was created
Definition at line 138 of file LinearSystem.hpp.
Referenced by LinearSystem(), Solve(), and ~LinearSystem().
|
private |
Stores one more than the highest index stored locally.
Definition at line 85 of file LinearSystem.hpp.
Referenced by GetOwnershipRange(), LinearSystem(), and SetPrecondMatrixIsDifferentFromLhs().
|
private |
For parallel code. Stores lowest index of vectors and lowest row of matrix stored locally.
Definition at line 84 of file LinearSystem.hpp.
Referenced by GetOwnershipRange(), LinearSystem(), and SetPrecondMatrixIsDifferentFromLhs().
|
private |
Pointer to vector containing a list of bath nodes
Definition at line 116 of file LinearSystem.hpp.
Referenced by SetPcType(), and Solve().
|
private |
Stores a pointer to a purpose-build preconditioner
Definition at line 109 of file LinearSystem.hpp.
Referenced by SetPcType(), Solve(), and ~LinearSystem().
|
private |
Context for KSPDefaultConverged()
Definition at line 135 of file LinearSystem.hpp.
Referenced by Solve(), and ~LinearSystem().
|
private |
Preconditioner type (see PETSc PCSetType() )
Definition at line 103 of file LinearSystem.hpp.
Referenced by LinearSystem(), serialize(), SetPcType(), and Solve().
|
private |
Stores a pointer to a purpose-build preconditioner
Definition at line 111 of file LinearSystem.hpp.
Referenced by SetPcType(), Solve(), and ~LinearSystem().
|
private |
The matrix used for preconditioning.
Definition at line 76 of file LinearSystem.hpp.
Referenced by FinalisePrecondMatrix(), rGetPrecondMatrix(), SetPrecondMatrixIsDifferentFromLhs(), Solve(), and ~LinearSystem().
|
private |
Whether the matrix used for preconditioning is the same as the LHS
Definition at line 119 of file LinearSystem.hpp.
Referenced by rGetPrecondMatrix(), SetPrecondMatrixIsDifferentFromLhs(), Solve(), and ~LinearSystem().
|
private |
Stores a pointer to a purpose-build preconditioner
Definition at line 113 of file LinearSystem.hpp.
Referenced by SetPcType(), Solve(), and ~LinearSystem().
|
private |
The right-hand side vector.
Definition at line 77 of file LinearSystem.hpp.
Referenced by AddRhsMultipleValues(), AddToRhsVectorElement(), DisplayRhs(), FinaliseRhsVector(), GetRhsVector(), GetRhsVectorElement(), LinearSystem(), rGetRhsVector(), SetRhsVectorElement(), Solve(), ZeroRhsVector(), and ~LinearSystem().
|
private |
The max number of nonzero entries expected on a LHS row
Definition at line 122 of file LinearSystem.hpp.
Referenced by LinearSystem(), and SetPrecondMatrixIsDifferentFromLhs().
|
private |
The size of the linear system.
Definition at line 78 of file LinearSystem.hpp.
Referenced by GetSize(), LinearSystem(), SetPrecondMatrixIsDifferentFromLhs(), and Solve().
|
private |
absolute or relative tolerance of the KSP solver
Definition at line 96 of file LinearSystem.hpp.
Referenced by serialize(), SetAbsoluteTolerance(), SetRelativeTolerance(), and Solve().
|
private |
Sets either absolute or relative tolerance of the KSP solver. Default is to false
Definition at line 101 of file LinearSystem.hpp.
Referenced by serialize(), SetAbsoluteTolerance(), SetRelativeTolerance(), and Solve().
|
private |
Whether to use fixed number of iterations
Definition at line 125 of file LinearSystem.hpp.
Referenced by SetUseFixedNumberIterations(), and Solve().