LinearSystem Class Reference

#include <LinearSystem.hpp>

Collaboration diagram for LinearSystem:
Collaboration graph
[legend]

List of all members.

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 ()
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::LinearSystem ( PetscInt  lhsVectorSize,
unsigned  rowPreallocation = UINT_MAX 
)

Constructor.

Parameters:
lhsVectorSize the size of the LHS vector
rowPreallocation the 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::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:
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

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::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:
residualVector the residual vector
jacobianMatrix the 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::LinearSystem ( PetscInt  lhsVectorSize,
Mat  lhsMatrix,
Vec  rhsVector 
)

Alternative constructor for archiving.

Parameters:
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 (  ) 

Member Function Documentation

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

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

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

Add the specified value to an entry of the matrix.

Parameters:
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.

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

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

Parameters:
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 (  ) 
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 (  ) 
void LinearSystem::FinalisePrecondMatrix (  ) 

Sets up the PETSc matrix used for preconditioning.

Definition at line 342 of file LinearSystem.cpp.

References PetscMatTools::Finalise(), and mPrecondMatrix.

Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem().

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

Definition at line 513 of file LinearSystem.cpp.

References mLhsMatrix.

Referenced by AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseForBath().

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

Definition at line 478 of file LinearSystem.cpp.

References PetscMatTools::GetElement(), and mLhsMatrix.

Referenced by AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseForBath().

Vec LinearSystem::GetMatrixRowDistributed ( unsigned  rowIndex  ) 

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

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

Definition at line 377 of file LinearSystem.cpp.

References mLhsMatrix.

Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().

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

Definition at line 488 of file LinearSystem.cpp.

References mKspIsSetup, and mKspSolver.

void LinearSystem::GetOwnershipRange ( PetscInt lo,
PetscInt hi 
)
Returns:
this process's ownership range of the contents of the system.
Parameters:
lo lowest index owned by this process
hi highest index owned by this process

Definition at line 472 of file LinearSystem.cpp.

References mOwnershipRangeHi, and mOwnershipRangeLo.

Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem(), and AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseForBath().

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

Definition at line 503 of file LinearSystem.cpp.

References mRhsVector.

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

Definition at line 483 of file LinearSystem.cpp.

References PetscVecTools::GetElement(), and mRhsVector.

unsigned LinearSystem::GetSize (  )  const
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 557 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, mMatNullSpace, and PETSC_DESTROY_PARAM.

Referenced by BidomainProblem< DIM >::AtBeginningOfTimestep().

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 1216 of file LinearSystem.cpp.

References mForceSpectrumReevaluation, mKspIsSetup, mKspSolver, and PETSC_DESTROY_PARAM.

Referenced by AbstractDynamicLinearPdeSolver< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve().

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

Should only be used by the BoundaryConditionsContainer.

Definition at line 528 of file LinearSystem.cpp.

References mDirichletBoundaryConditionsVector.

Referenced by BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem().

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

Definition at line 518 of file LinearSystem.cpp.

References EXCEPTION, mPrecondMatrix, and mPrecondMatrixIsNotLhs.

Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem().

Vec & LinearSystem::rGetRhsVector (  ) 
template<class Archive >
void LinearSystem::serialize ( Archive &  archive,
const unsigned int  version 
) [inline, private]

Archive the member variables.

Parameters:
archive 
version 

Definition at line 162 of file LinearSystem.hpp.

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

void LinearSystem::SetAbsoluteTolerance ( double  absoluteTolerance  ) 
void LinearSystem::SetKspType ( const char *  kspType  ) 
void LinearSystem::SetMatrixElement ( PetscInt  row,
PetscInt  col,
double  value 
)

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

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

Referenced by AbstractExtendedBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseLinearSystem(), and AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseLinearSystem().

void LinearSystem::SetMatrixIsConstant ( bool  matrixIsConstant  ) 

Set mMatrixIsConstant.

Parameters:
matrixIsConstant whether the matrix is constant

Definition at line 568 of file LinearSystem.cpp.

References mMatrixIsConstant.

void LinearSystem::SetMatrixIsSymmetric ( bool  isSymmetric = true  ) 
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:
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().

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

Referenced by AbstractExtendedBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseLinearSystem(), and AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseLinearSystem().

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:
pcType the preconditioner type
pBathNodes the 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 603 of file LinearSystem.cpp.

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

Referenced by OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), AbstractExtendedBidomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), and AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve().

void LinearSystem::SetPrecondMatrixIsDifferentFromLhs ( bool  precondIsDifferent = true  ) 

Set method for mPrecondMatrixIsNotLhs.

Parameters:
precondIsDifferent whether the matrix used for preconditioning is the same as the LHS.

Definition at line 1184 of file LinearSystem.cpp.

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

Referenced by MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::SetupLinearSystem().

void LinearSystem::SetRelativeTolerance ( double  relativeTolerance  ) 
void LinearSystem::SetRhsVectorElement ( PetscInt  row,
double  value 
)
void LinearSystem::SetUseFixedNumberIterations ( bool  useFixedNumberIterations = true,
unsigned  evaluateNumItsEveryNSolves = UINT_MAX 
)

Set method for mUseFixedNumberIterations

Parameters:
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 1209 of file LinearSystem.cpp.

References mEvaluateNumItsEveryNSolves, and mUseFixedNumberIterations.

Referenced by OperatorSplittingMonodomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), MonodomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), and AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve().

Vec LinearSystem::Solve ( Vec  lhsGuess = NULL  ) 
void LinearSystem::SwitchWriteModeLhsMatrix (  ) 
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.

Parameters:
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.

Parameters:
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

Parameters:
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(), BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyPeriodicBcsToLinearProblem(), AbstractExtendedBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseLinearSystem(), and AbstractBidomainSolver< ELEMENT_DIM, SPACE_DIM >::FinaliseLinearSystem().

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 Function Documentation

friend class boost::serialization::access [friend]

Needed for serialization.

Definition at line 154 of file LinearSystem.hpp.


Member Data Documentation

Whether we need to destroy the PETSc matrix and vector in our destructor

Definition at line 90 of file LinearSystem.hpp.

Referenced by ~LinearSystem().

Storage for efficient application of Dirichlet BCs, see AbstractBoundaryConditionsContainer

Definition at line 105 of file LinearSystem.hpp.

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

PetscReal LinearSystem::mEigMax [private]

Preconditioned operator largest eigenvalue

Definition at line 144 of file LinearSystem.hpp.

Referenced by Solve().

PetscReal LinearSystem::mEigMin [private]

Preconditioned operator smallest eigenvalue

Definition at line 141 of file LinearSystem.hpp.

Referenced by Solve().

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

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

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

KSP LinearSystem::mKspSolver [private]
std::string LinearSystem::mKspType [private]

KSP solver type (see PETSc KSPSetType() )

Definition at line 102 of file LinearSystem.hpp.

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

MatNullSpace LinearSystem::mMatNullSpace [private]

PETSc null matrix.

Definition at line 87 of file LinearSystem.hpp.

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

Whether the matrix is unchanged each time Solve() is called

Definition at line 95 of file LinearSystem.hpp.

Referenced by serialize(), SetMatrixIsConstant(), and Solve().

Yes, it really is stored as a double.

Definition at line 94 of file LinearSystem.hpp.

Referenced by serialize(), and Solve().

Number of solves performed since the current object was created

Definition at line 138 of file LinearSystem.hpp.

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

Stores one more than the highest index stored locally.

Definition at line 85 of file LinearSystem.hpp.

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

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

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

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 SetPcType(), Solve(), and ~LinearSystem().

Context for KSPDefaultConverged()

Definition at line 135 of file LinearSystem.hpp.

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

std::string LinearSystem::mPcType [private]

Preconditioner type (see PETSc PCSetType() )

Definition at line 103 of file LinearSystem.hpp.

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

Stores a pointer to a purpose-build preconditioner

Definition at line 111 of file LinearSystem.hpp.

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

The matrix used for preconditioning.

Definition at line 76 of file LinearSystem.hpp.

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

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

Stores a pointer to a purpose-build preconditioner

Definition at line 113 of file LinearSystem.hpp.

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

The max number of nonzero entries expected on a LHS row

Definition at line 122 of file LinearSystem.hpp.

Referenced by LinearSystem(), and SetPrecondMatrixIsDifferentFromLhs().

The size of the linear system.

Definition at line 78 of file LinearSystem.hpp.

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

absolute or relative tolerance of the KSP solver

Definition at line 96 of file LinearSystem.hpp.

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

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

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:

Generated by  doxygen 1.6.2