Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
LinearSystem Class Reference

#include <LinearSystem.hpp>

+ Collaboration diagram for LinearSystem:

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=nullptr)
 
void SetRhsVectorElement (PetscInt row, double value)
 
void AddToRhsVectorElement (PetscInt row, double value)
 
unsigned GetSize () const
 
void SetNullBasis (Vec nullbasis[], unsigned numberOfBases)
 
void RemoveNullSpace ()
 
VecrGetRhsVector ()
 
Vec GetRhsVector () const
 
MatrGetLhsMatrix ()
 
Mat GetLhsMatrix () const
 
MatrGetPrecondMatrix ()
 
VecrGetDirichletBoundaryConditionsVector ()
 
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
 
PCBlockDiagonalmpBlockDiagonalPC
 
PCLDUFactorisationmpLDUFactorisationPC
 
PCTwoLevelsBlockDiagonalmpTwoLevelsBlockDiagonalPC
 
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
 

Detailed Description

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 & Destructor Documentation

◆ LinearSystem() [1/4]

LinearSystem::LinearSystem ( PetscInt  lhsVectorSize,
unsigned  rowPreallocation = UINT_MAX 
)

Constructor.

Parameters
lhsVectorSizethe size of the LHS vector
rowPreallocationthe max number of nonzero entries expected on a row
  • a value of 0 is allowed: no preallocation is then done and the user must preallocate the memory for the matrix themselves.
  • the default value allows for small size systems to be set as dense matrices automatically.
Todo:
: if we create a linear system object outside a cardiac solver, these are gonna be the default solver and preconditioner. Not consitent with ChasteDefaults.xml though...

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() [2/4]

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.

Parameters
templateVectora PETSc vec
rowPreallocationthe max number of nonzero entries expected on a row
newAllocationErrortells 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
Todo:
: if we create a linear system object outside a cardiac solver, these are gonna be the default solver and preconditioner. Not consistent with ChasteDefaults.xml though...

Definition at line 146 of file LinearSystem.cpp.

References mKspType, mLhsMatrix, mNumSolves, mOwnershipRangeHi, mOwnershipRangeLo, mPcType, mRhsVector, mRowPreallocation, mSize, and PetscTools::SetupMat().

◆ LinearSystem() [3/4]

LinearSystem::LinearSystem ( Vec  residualVector,
Mat  jacobianMatrix 
)

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.

Parameters
residualVectorthe residual vector
jacobianMatrixthe Jacobian matrix
Todo:
: if we create a linear system object outside a cardiac solver, these are gonna be the default solver and preconditioner. Not consitent with ChasteDefaults.xml though...

Definition at line 187 of file LinearSystem.cpp.

References mKspType, mLhsMatrix, mNumSolves, mOwnershipRangeHi, mOwnershipRangeLo, mPcType, mRhsVector, mRowPreallocation, and mSize.

◆ LinearSystem() [4/4]

LinearSystem::LinearSystem ( PetscInt  lhsVectorSize,
Mat  lhsMatrix,
Vec  rhsVector 
)

Alternative constructor for archiving.

Parameters
lhsVectorSizethe size of the LHS vector
lhsMatrixthe RHS matrix
rhsVectorthe RHS vector

Definition at line 109 of file LinearSystem.cpp.

References mLhsMatrix, mNumSolves, mOwnershipRangeHi, mOwnershipRangeLo, and mRhsVector.

◆ ~LinearSystem()

Member Function Documentation

◆ AddLhsMultipleValues()

template<size_t MATRIX_SIZE>
void LinearSystem::AddLhsMultipleValues ( unsigned matrixRowAndColIndices,
c_matrix< double, MATRIX_SIZE, MATRIX_SIZE > &  rSmallMatrix 
)
inline

Add multiple values to the matrix of linear system.

Parameters
matrixRowAndColIndicesmapping from index of the ublas matrix (see param below) to index of the PETSc matrix of this linear system
rSmallMatrixUblas 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.

◆ AddRhsMultipleValues()

template<size_t VECTOR_SIZE>
void LinearSystem::AddRhsMultipleValues ( unsigned vectorIndices,
c_vector< double, VECTOR_SIZE > &  smallVector 
)
inline

Add multiple values to the RHS vector.

Parameters
vectorIndicesmapping from index of the ublas vector (see param below) to index of the vector of this linear system
smallVectorUblas 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.

◆ AddToMatrixElement()

void LinearSystem::AddToMatrixElement ( PetscInt  row,
PetscInt  col,
double  value 
)

Add the specified value to an entry of the matrix.

Parameters
rowthe row index
colthe column index
valuethe value for this entry

Definition at line 315 of file LinearSystem.cpp.

References PetscMatTools::AddToElement(), and mLhsMatrix.

◆ AddToRhsVectorElement()

void LinearSystem::AddToRhsVectorElement ( PetscInt  row,
double  value 
)

Add a value to an element of the right-hand side vector.

Parameters
rowthe row index
valuethe value to set this entry

Definition at line 357 of file LinearSystem.cpp.

References PetscVecTools::AddToElement(), and mRhsVector.

◆ AssembleFinalLinearSystem()

void LinearSystem::AssembleFinalLinearSystem ( )

Call this before Solve().

This calls FinaliseLhsMatrix() and FinaliseRhsVector().

Definition at line 320 of file LinearSystem.cpp.

References FinaliseLhsMatrix(), and FinaliseRhsVector().

◆ AssembleIntermediateLinearSystem()

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().

◆ DisplayMatrix()

void LinearSystem::DisplayMatrix ( )

Display the left-hand side matrix.

Definition at line 362 of file LinearSystem.cpp.

References PetscMatTools::Display(), and mLhsMatrix.

◆ DisplayRhs()

void LinearSystem::DisplayRhs ( )

Display the right-hand side vector.

Definition at line 367 of file LinearSystem.cpp.

References PetscVecTools::Display(), and mRhsVector.

◆ FinaliseLhsMatrix()

void LinearSystem::FinaliseLhsMatrix ( )

◆ FinalisePrecondMatrix()

void LinearSystem::FinalisePrecondMatrix ( )

Sets up the PETSc matrix used for preconditioning.

Definition at line 342 of file LinearSystem.cpp.

References PetscMatTools::Finalise(), and mPrecondMatrix.

◆ FinaliseRhsVector()

void LinearSystem::FinaliseRhsVector ( )

◆ GetLhsMatrix()

Mat LinearSystem::GetLhsMatrix ( ) const
Returns
access to the LHS matrix for archiving

Definition at line 524 of file LinearSystem.cpp.

References mLhsMatrix.

◆ GetMatrixElement()

double LinearSystem::GetMatrixElement ( PetscInt  row,
PetscInt  col 
)
Returns
an element of the matrix. May only be called for elements you own.
Parameters
rowthe row index
colthe column index

Definition at line 489 of file LinearSystem.cpp.

References PetscMatTools::GetElement(), and mLhsMatrix.

◆ GetMatrixRowDistributed()

Vec LinearSystem::GetMatrixRowDistributed ( unsigned  rowIndex)

Returns the i-th row of the LHS matrix as a distributed PETSc Vec

Parameters
rowIndexthe row index
Returns
rowIndex-th row of the matrix in distributed format

Definition at line 377 of file LinearSystem.cpp.

References PetscMatTools::GetMatrixRowDistributed(), and mLhsMatrix.

◆ GetNumIterations()

unsigned LinearSystem::GetNumIterations ( ) const
Returns
the number of iterations taken by the last Solve()

Definition at line 499 of file LinearSystem.cpp.

References mKspIsSetup, and mKspSolver.

◆ GetOwnershipRange()

void LinearSystem::GetOwnershipRange ( PetscInt lo,
PetscInt hi 
)

Return this process's ownership range of the contents of the system via arguments.

Parameters
lolowest index owned by this process (returned)
hihighest index owned by this process (returned)

Definition at line 483 of file LinearSystem.cpp.

References mOwnershipRangeHi, and mOwnershipRangeLo.

◆ GetRhsVector()

Vec LinearSystem::GetRhsVector ( ) const
Returns
access to the RHS vector for archiving

Definition at line 514 of file LinearSystem.cpp.

References mRhsVector.

◆ GetRhsVectorElement()

double LinearSystem::GetRhsVectorElement ( PetscInt  row)
Returns
an element of the RHS vector. May only be called for elements you own.
Parameters
rowthe row index

Definition at line 494 of file LinearSystem.cpp.

References PetscVecTools::GetElement(), and mRhsVector.

◆ GetSize()

unsigned LinearSystem::GetSize ( ) const
Returns
mSize.

Definition at line 413 of file LinearSystem.cpp.

References mSize.

◆ IsMatrixSymmetric()

bool LinearSystem::IsMatrixSymmetric ( )

Get whether PETSc considers the matrix in this linear system as symmetric or not.

Returns
whether the matrix is symmetric or not.

Definition at line 568 of file LinearSystem.cpp.

References mLhsMatrix.

◆ RemoveNullSpace()

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().

◆ 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).

Todo:
#1695 Store this number in a member variable.

Definition at line 1249 of file LinearSystem.cpp.

References mForceSpectrumReevaluation, mKspIsSetup, mKspSolver, PETSC_DESTROY_PARAM, and PetscTools::SetOption().

Referenced by RemoveNullSpace().

◆ rGetDirichletBoundaryConditionsVector()

Vec & LinearSystem::rGetDirichletBoundaryConditionsVector ( )
Returns
access to the dirichlet boundary conditions vector.

Should only be used by the BoundaryConditionsContainer.

Definition at line 539 of file LinearSystem.cpp.

References mDirichletBoundaryConditionsVector.

◆ rGetLhsMatrix()

◆ rGetPrecondMatrix()

Mat & LinearSystem::rGetPrecondMatrix ( )
Returns
access to the matrix used for preconditioning.

Definition at line 529 of file LinearSystem.cpp.

References EXCEPTION, mPrecondMatrix, and mPrecondMatrixIsNotLhs.

◆ rGetRhsVector()

◆ serialize()

template<class Archive >
void LinearSystem::serialize ( Archive &  archive,
const unsigned int  version 
)
inlineprivate

Archive the member variables.

Parameters
archive
version

Definition at line 162 of file LinearSystem.hpp.

References mKspType, mMatrixIsConstant, mNonZerosUsed, mPcType, mTolerance, and mUseAbsoluteTolerance.

◆ SetAbsoluteTolerance()

void LinearSystem::SetAbsoluteTolerance ( double  absoluteTolerance)

Set the absolute tolerance.

Parameters
absoluteTolerancethe absolute tolerance

Definition at line 594 of file LinearSystem.cpp.

References mKspIsSetup, mKspSolver, mTolerance, and mUseAbsoluteTolerance.

◆ SetKspType()

void LinearSystem::SetKspType ( const char *  kspType)

Set the KSP solver type (see PETSc KSPSetType() for valid arguments).

Parameters
kspTypethe KSP solver type

Definition at line 604 of file LinearSystem.cpp.

References mKspIsSetup, mKspSolver, and mKspType.

◆ SetMatrixElement()

void LinearSystem::SetMatrixElement ( PetscInt  row,
PetscInt  col,
double  value 
)

Change one of the entires of the matrix to the specified value.

Parameters
rowthe row index
colthe column index
valuethe value for this entry

Definition at line 310 of file LinearSystem.cpp.

References mLhsMatrix, and PetscMatTools::SetElement().

◆ SetMatrixIsConstant()

void LinearSystem::SetMatrixIsConstant ( bool  matrixIsConstant)

Set mMatrixIsConstant.

Parameters
matrixIsConstantwhether the matrix is constant

Definition at line 579 of file LinearSystem.cpp.

References mMatrixIsConstant.

◆ SetMatrixIsSymmetric()

void LinearSystem::SetMatrixIsSymmetric ( bool  isSymmetric = true)

Force PETSc to treat the matrix in this linear system as symmetric from now on.

Parameters
isSymmetricwhether the matrix is symmetric or not (defaults to true)
Todo:
: shall we allow modifying the symmetry flag anytime?

Definition at line 544 of file LinearSystem.cpp.

References mLhsMatrix, and PetscMatTools::SetOption().

◆ SetMatrixRow()

void LinearSystem::SetMatrixRow ( PetscInt  row,
double  value 
)

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)

Parameters
rowthe row index
valuethe value to set each entry in this row

Definition at line 372 of file LinearSystem.cpp.

References mLhsMatrix, and PetscMatTools::SetRow().

◆ SetNullBasis()

void LinearSystem::SetNullBasis ( Vec  nullbasis[],
unsigned  numberOfBases 
)

Set the null basis of the linear system. In debug mode we test for orthonormality and throw EXCEPTIONs:

  • each direction is of unit length (to within a tolerance)
  • the directions are mutually orthogonal (to within a tolerance) Note that in NDEBUG (optimised mode) these tests are skipped and the parameters are passed directly into the LinearSystem.
Parameters
nullbasisan array PETSc vectors containing orthogonal directions in the nullspace
numberOfBasesthe number of directions (size of nullbasis array)

Definition at line 418 of file LinearSystem.cpp.

References EXCEPTION, and mMatNullSpace.

◆ SetPcType()

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).

Parameters
pcTypethe preconditioner type
pBathNodesthe list of nodes defining the bath
Todo:
: #1082 is this the way of defining a null pointer as the default value of pBathNodes?
Todo:
: #1082 use a single pointer to abstract class
Todo:
: #1082 use a single pointer to abstract class
Todo:
: #1082 use a single pointer to abstract class

Definition at line 614 of file LinearSystem.cpp.

References mKspIsSetup, mKspSolver, mpBathNodes, mpBlockDiagonalPC, mPcType, mpLDUFactorisationPC, mpTwoLevelsBlockDiagonalPC, and TERMINATE.

◆ SetPrecondMatrixIsDifferentFromLhs()

void LinearSystem::SetPrecondMatrixIsDifferentFromLhs ( bool  precondIsDifferent = true)

Set method for mPrecondMatrixIsNotLhs.

Parameters
precondIsDifferentwhether the matrix used for preconditioning is the same as the LHS.

Definition at line 1217 of file LinearSystem.cpp.

References mOwnershipRangeHi, mOwnershipRangeLo, mPrecondMatrix, mPrecondMatrixIsNotLhs, mRowPreallocation, mSize, NEVER_REACHED, and PetscTools::SetupMat().

◆ SetRelativeTolerance()

void LinearSystem::SetRelativeTolerance ( double  relativeTolerance)

Set the relative tolerance.

Parameters
relativeTolerancethe relative tolerance

Definition at line 584 of file LinearSystem.cpp.

References mKspIsSetup, mKspSolver, mTolerance, and mUseAbsoluteTolerance.

◆ SetRhsVectorElement()

void LinearSystem::SetRhsVectorElement ( PetscInt  row,
double  value 
)

Set an element of the right-hand side vector to a given value.

Parameters
rowthe row index
valuethe 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 >::ApplyPeriodicBcsToLinearProblem().

◆ SetUseFixedNumberIterations()

void LinearSystem::SetUseFixedNumberIterations ( bool  useFixedNumberIterations = true,
unsigned  evaluateNumItsEveryNSolves = UINT_MAX 
)

Set method for mUseFixedNumberIterations

Parameters
useFixedNumberIterationswhether to use fixed number of iterations
evaluateNumItsEveryNSolvestells 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 1242 of file LinearSystem.cpp.

References mEvaluateNumItsEveryNSolves, and mUseFixedNumberIterations.

◆ Solve()

Vec LinearSystem::Solve ( Vec  lhsGuess = nullptr)

◆ SwitchWriteModeLhsMatrix()

void LinearSystem::SwitchWriteModeLhsMatrix ( )

◆ ZeroLhsMatrix()

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().

◆ 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().

◆ ZeroMatrixColumn()

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.

Parameters
colthe column index

Definition at line 392 of file LinearSystem.cpp.

References mLhsMatrix, and PetscMatTools::ZeroColumn().

◆ ZeroMatrixRowsAndColumnsWithValueOnDiagonal()

void LinearSystem::ZeroMatrixRowsAndColumnsWithValueOnDiagonal ( std::vector< unsigned > &  rRowColIndices,
double  diagonalValue 
)

Zero several rows and columns of the matrix, putting a given value on the diagonal.

Parameters
rRowColIndicesA list of indices. All the rows with these indices, and all the columns with these indices, will be zeroed.
diagonalValuevalue to put in the diagonal entries (of the zeroed rows)

Definition at line 387 of file LinearSystem.cpp.

References mLhsMatrix, and PetscMatTools::ZeroRowsAndColumnsWithValueOnDiagonal().

◆ ZeroMatrixRowsWithValueOnDiagonal()

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

Parameters
rRowsstd::vector of rows to be zeroed
diagonalValuevalue 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 >::ApplyPeriodicBcsToLinearProblem().

◆ ZeroRhsVector()

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().

Friends And Related Symbol Documentation

◆ boost::serialization::access

friend class boost::serialization::access
friend

Needed for serialization.

Definition at line 154 of file LinearSystem.hpp.

◆ TestChebyshevIteration

friend class TestChebyshevIteration
friend

Definition at line 71 of file LinearSystem.hpp.

◆ TestLinearSystem

friend class TestLinearSystem
friend

Definition at line 67 of file LinearSystem.hpp.

◆ TestPCBlockDiagonal

friend class TestPCBlockDiagonal
friend

Definition at line 68 of file LinearSystem.hpp.

◆ TestPCLDUFactorisation

friend class TestPCLDUFactorisation
friend

Definition at line 70 of file LinearSystem.hpp.

◆ TestPCTwoLevelsBlockDiagonal

friend class TestPCTwoLevelsBlockDiagonal
friend

Definition at line 69 of file LinearSystem.hpp.

Member Data Documentation

◆ mDestroyMatAndVec

bool LinearSystem::mDestroyMatAndVec
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().

◆ mDirichletBoundaryConditionsVector

Vec LinearSystem::mDirichletBoundaryConditionsVector
private

Storage for efficient application of Dirichlet BCs, see AbstractBoundaryConditionsContainer

Definition at line 105 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), and rGetDirichletBoundaryConditionsVector().

◆ mEigMax

PetscReal LinearSystem::mEigMax
private

Preconditioned operator largest eigenvalue

Definition at line 144 of file LinearSystem.hpp.

Referenced by Solve().

◆ mEigMin

PetscReal LinearSystem::mEigMin
private

Preconditioned operator smallest eigenvalue

Definition at line 141 of file LinearSystem.hpp.

Referenced by Solve().

◆ mEvaluateNumItsEveryNSolves

unsigned LinearSystem::mEvaluateNumItsEveryNSolves
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().

◆ mForceSpectrumReevaluation

bool LinearSystem::mForceSpectrumReevaluation
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().

◆ mKspIsSetup

bool LinearSystem::mKspIsSetup
private

Used by Solve method to track whether KSP has been used.

Definition at line 93 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), GetNumIterations(), RemoveNullSpace(), ResetKspSolver(), SetAbsoluteTolerance(), SetKspType(), SetPcType(), SetRelativeTolerance(), and Solve().

◆ mKspSolver

KSP LinearSystem::mKspSolver
private

◆ mKspType

std::string LinearSystem::mKspType
private

KSP solver type (see PETSc KSPSetType() )

Definition at line 102 of file LinearSystem.hpp.

Referenced by LinearSystem(), LinearSystem(), LinearSystem(), serialize(), SetKspType(), and Solve().

◆ mLhsMatrix

◆ mMatNullSpace

MatNullSpace LinearSystem::mMatNullSpace
private

PETSc null matrix.

Definition at line 87 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), RemoveNullSpace(), SetNullBasis(), and Solve().

◆ mMatrixIsConstant

bool LinearSystem::mMatrixIsConstant
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().

◆ mNonZerosUsed

double LinearSystem::mNonZerosUsed
private

Yes, it really is stored as a double.

Definition at line 94 of file LinearSystem.hpp.

Referenced by serialize(), and Solve().

◆ mNumSolves

unsigned LinearSystem::mNumSolves
private

Number of solves performed since the current object was created

Definition at line 138 of file LinearSystem.hpp.

Referenced by LinearSystem(), LinearSystem(), LinearSystem(), LinearSystem(), ~LinearSystem(), and Solve().

◆ mOwnershipRangeHi

PetscInt LinearSystem::mOwnershipRangeHi
private

Stores one more than the highest index stored locally.

Definition at line 85 of file LinearSystem.hpp.

Referenced by LinearSystem(), LinearSystem(), LinearSystem(), LinearSystem(), GetOwnershipRange(), and SetPrecondMatrixIsDifferentFromLhs().

◆ mOwnershipRangeLo

PetscInt LinearSystem::mOwnershipRangeLo
private
Todo:
Verify claim that ownership range for Vec and Mat is same. This should only matter for efficiency if the claim is false.

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 LinearSystem(), LinearSystem(), LinearSystem(), LinearSystem(), GetOwnershipRange(), and SetPrecondMatrixIsDifferentFromLhs().

◆ mpBathNodes

boost::shared_ptr<std::vector<PetscInt> > LinearSystem::mpBathNodes
private

Pointer to vector containing a list of bath nodes

Definition at line 116 of file LinearSystem.hpp.

Referenced by SetPcType(), and Solve().

◆ mpBlockDiagonalPC

PCBlockDiagonal* LinearSystem::mpBlockDiagonalPC
private
Todo:
: #1082 Create an abstract class for the preconditioners and use a single pointer

Stores a pointer to a purpose-build preconditioner

Definition at line 109 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), SetPcType(), and Solve().

◆ mpConvergenceTestContext

void* LinearSystem::mpConvergenceTestContext
private

Context for KSPDefaultConverged()

Definition at line 135 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), and Solve().

◆ mPcType

std::string LinearSystem::mPcType
private

Preconditioner type (see PETSc PCSetType() )

Definition at line 103 of file LinearSystem.hpp.

Referenced by LinearSystem(), LinearSystem(), LinearSystem(), serialize(), SetPcType(), and Solve().

◆ mpLDUFactorisationPC

PCLDUFactorisation* LinearSystem::mpLDUFactorisationPC
private

Stores a pointer to a purpose-build preconditioner

Definition at line 111 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), SetPcType(), and Solve().

◆ mPrecondMatrix

Mat LinearSystem::mPrecondMatrix
private

The matrix used for preconditioning.

Definition at line 76 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), FinalisePrecondMatrix(), rGetPrecondMatrix(), SetPrecondMatrixIsDifferentFromLhs(), and Solve().

◆ mPrecondMatrixIsNotLhs

bool LinearSystem::mPrecondMatrixIsNotLhs
private

Whether the matrix used for preconditioning is the same as the LHS

Definition at line 119 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), rGetPrecondMatrix(), SetPrecondMatrixIsDifferentFromLhs(), and Solve().

◆ mpTwoLevelsBlockDiagonalPC

PCTwoLevelsBlockDiagonal* LinearSystem::mpTwoLevelsBlockDiagonalPC
private

Stores a pointer to a purpose-build preconditioner

Definition at line 113 of file LinearSystem.hpp.

Referenced by ~LinearSystem(), SetPcType(), and Solve().

◆ mRhsVector

◆ mRowPreallocation

unsigned LinearSystem::mRowPreallocation
private

The max number of nonzero entries expected on a LHS row

Definition at line 122 of file LinearSystem.hpp.

Referenced by LinearSystem(), LinearSystem(), LinearSystem(), and SetPrecondMatrixIsDifferentFromLhs().

◆ mSize

PetscInt LinearSystem::mSize
private

The size of the linear system.

Definition at line 78 of file LinearSystem.hpp.

Referenced by LinearSystem(), LinearSystem(), LinearSystem(), GetSize(), SetPrecondMatrixIsDifferentFromLhs(), and Solve().

◆ mTolerance

double LinearSystem::mTolerance
private

absolute or relative tolerance of the KSP solver

Definition at line 96 of file LinearSystem.hpp.

Referenced by serialize(), SetAbsoluteTolerance(), SetRelativeTolerance(), and Solve().

◆ mUseAbsoluteTolerance

bool LinearSystem::mUseAbsoluteTolerance
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().

◆ mUseFixedNumberIterations

bool LinearSystem::mUseFixedNumberIterations
private

Whether to use fixed number of iterations

Definition at line 125 of file LinearSystem.hpp.

Referenced by SetUseFixedNumberIterations(), and Solve().


The documentation for this class was generated from the following files: