LinearSystem Class Reference

#include <LinearSystem.hpp>

Collaboration diagram for LinearSystem:

Collaboration graph
[legend]

List of all members.

Public Member Functions

 LinearSystem (PetscInt lhsVectorSize, MatType matType=(MatType) MATMPIAIJ)
 LinearSystem (Vec templateVector)
 LinearSystem (Vec residualVector, Mat jacobianMatrix)
 LinearSystem (PetscInt lhsVectorSize, Mat lhsMatrix, Vec rhsVector, MatType matType=(MatType) MATMPIAIJ)
 ~LinearSystem ()
void SetupVectorAndMatrix (MatType matType)
void SetMatrixElement (PetscInt row, PetscInt col, double value)
void AddToMatrixElement (PetscInt row, PetscInt col, double value)
void AssembleFinalLinearSystem ()
void AssembleIntermediateLinearSystem ()
void AssembleFinalLhsMatrix ()
void AssembleIntermediateLhsMatrix ()
void AssembleRhsVector ()
void SetMatrixIsSymmetric (bool isSymmetric=true)
void SetMatrixIsConstant (bool matrixIsConstant)
void SetRelativeTolerance (double relativeTolerance)
void SetAbsoluteTolerance (double absoluteTolerance)
void SetKspType (const char *kspType)
void SetPcType (const char *pcType)
void DisplayMatrix ()
void DisplayRhs ()
void SetMatrixRow (PetscInt row, double value)
void ZeroMatrixRow (PetscInt row)
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)
Vec & rGetRhsVector ()
Vec GetRhsVector () const
Mat & rGetLhsMatrix ()
Mat GetLhsMatrix () const
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)

Private Member Functions

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

Private Attributes

Mat mLhsMatrix
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

Friends

class TestLinearSystem
class TestPCBlockDiagonal
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 57 of file LinearSystem.hpp.


Constructor & Destructor Documentation

LinearSystem::LinearSystem ( PetscInt  lhsVectorSize,
MatType  matType = (MatType) MATMPIAIJ 
)

Constructor.

Parameters:
lhsVectorSize the size of the LHS vector
matType the type of matrix (defaults to MATMPIAIJ)

Todo:
: if we create a linear system object outside a cardiac assembler, these are gonna be the default solver and preconditioner. Not consitent with ChasteDefaults.xml though...

Definition at line 42 of file LinearSystem.cpp.

References mKspType, mPcType, and SetupVectorAndMatrix().

LinearSystem::LinearSystem ( Vec  templateVector  ) 

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

Todo:
: if we create a linear system object outside a cardiac assembler, these are gonna be the default solver and preconditioner. Not consitent with ChasteDefaults.xml though...

Definition at line 93 of file LinearSystem.cpp.

References mKspType, mLhsMatrix, mOwnershipRangeHi, mOwnershipRangeLo, mPcType, mRhsVector, 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 assembler, these are gonna be the default solver and preconditioner. Not consitent with ChasteDefaults.xml though...

Definition at line 122 of file LinearSystem.cpp.

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

LinearSystem::LinearSystem ( PetscInt  lhsVectorSize,
Mat  lhsMatrix,
Vec  rhsVector,
MatType  matType = (MatType) MATMPIAIJ 
)

Alternative constructor for archiving.

Parameters:
lhsVectorSize the size of the LHS vector
lhsMatrix the RHS matrix
rhsVector the RHS vector
matType the type of matrix (defaults to MATMPIAIJ)

Definition at line 68 of file LinearSystem.cpp.

References mLhsMatrix, mOwnershipRangeHi, mOwnershipRangeLo, and mRhsVector.

LinearSystem::~LinearSystem (  ) 

Destructor.

Todo:
Never tested in linalg component

Definition at line 165 of file LinearSystem.cpp.

References mDestroyMatAndVec, mDirichletBoundaryConditionsVector, mKspIsSetup, mKspSolver, mLhsMatrix, mMatNullSpace, mpBlockDiagonalPC, and mRhsVector.


Member Function Documentation

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

Archive the member variables.

Parameters:
archive 
version 

Definition at line 112 of file LinearSystem.hpp.

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

void LinearSystem::SetupVectorAndMatrix ( MatType  matType  ) 

Helper method for the constructor. Initializes the LHS matrix and RHS vector.

Parameters:
matType the type of matrix

Definition at line 207 of file LinearSystem.cpp.

References mLhsMatrix, mOwnershipRangeHi, mOwnershipRangeLo, mRhsVector, mSize, and PetscTools::SetupMat().

Referenced by LinearSystem().

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

References mLhsMatrix, mOwnershipRangeHi, and mOwnershipRangeLo.

Referenced by AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem(), BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::FinaliseAssembleSystem(), and SetMatrixRow().

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

References mLhsMatrix, mOwnershipRangeHi, and mOwnershipRangeLo.

void LinearSystem::AssembleFinalLinearSystem (  ) 

void LinearSystem::AssembleIntermediateLinearSystem (  ) 

Should be called before AddToMatrixElement.

This calls AssembleIntermediateLhsMatrix() and AssembleRhsVector().

Definition at line 240 of file LinearSystem.cpp.

References AssembleIntermediateLhsMatrix(), and AssembleRhsVector().

void LinearSystem::AssembleFinalLhsMatrix (  ) 

void LinearSystem::AssembleIntermediateLhsMatrix (  ) 

void LinearSystem::AssembleRhsVector (  ) 

void LinearSystem::SetMatrixIsSymmetric ( bool  isSymmetric = true  ) 

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

Parameters:
isSymmetric whether the matrix is symmetric or not (defaults to true)

Definition at line 465 of file LinearSystem.cpp.

References mLhsMatrix.

Referenced by MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), and BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve().

void LinearSystem::SetMatrixIsConstant ( bool  matrixIsConstant  ) 

Set mMatrixIsConstant.

Parameters:
matrixIsConstant whether the matrix is constant

Definition at line 480 of file LinearSystem.cpp.

References mMatrixIsConstant.

Referenced by AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::InitialiseForSolve().

void LinearSystem::SetRelativeTolerance ( double  relativeTolerance  ) 

Set the relative tolerance.

Parameters:
relativeTolerance the relative tolerance

Definition at line 485 of file LinearSystem.cpp.

References mKspIsSetup, mKspSolver, mTolerance, and mUseAbsoluteTolerance.

Referenced by MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), and BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve().

void LinearSystem::SetAbsoluteTolerance ( double  absoluteTolerance  ) 

Set the absolute tolerance.

Parameters:
absoluteTolerance the absolute tolerance

Definition at line 495 of file LinearSystem.cpp.

References mKspIsSetup, mKspSolver, mTolerance, and mUseAbsoluteTolerance.

Referenced by MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), and BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve().

void LinearSystem::SetKspType ( const char *  kspType  ) 

void LinearSystem::SetPcType ( const char *  pcType  ) 

Set the preconditioner type (see PETSc PCSetType() for valid arguments).

Parameters:
pcType the preconditioner type

Definition at line 515 of file LinearSystem.cpp.

References mKspIsSetup, mKspSolver, mpBlockDiagonalPC, and mPcType.

Referenced by MonodomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve(), and BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::InitialiseForSolve().

void LinearSystem::DisplayMatrix (  ) 

Display the left-hand side matrix.

Definition at line 280 of file LinearSystem.cpp.

References mLhsMatrix.

void LinearSystem::DisplayRhs (  ) 

Display the right-hand side vector.

Definition at line 285 of file LinearSystem.cpp.

References mRhsVector.

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

References mLhsMatrix, mOwnershipRangeHi, mOwnershipRangeLo, and SetMatrixElement().

void LinearSystem::ZeroMatrixRow ( PetscInt  row  ) 

Zero a row of the left-hand side matrix. This method is a collective call (all processes should call it together). If processes call it with different arguments then its results may not be predictable.

Parameters:
row the row index

Definition at line 303 of file LinearSystem.cpp.

References mLhsMatrix.

Referenced by AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem(), and BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::FinaliseAssembleSystem().

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

References GetMatrixElement(), mLhsMatrix, mOwnershipRangeHi, and mOwnershipRangeLo.

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

void LinearSystem::ZeroLhsMatrix (  ) 

void LinearSystem::ZeroRhsVector (  ) 

void LinearSystem::ZeroLinearSystem (  ) 

Zero all entries of the left-hand side matrix and right-hand side vector.

Definition at line 379 of file LinearSystem.cpp.

References ZeroLhsMatrix(), and ZeroRhsVector().

Referenced by SimpleNewtonNonlinearSolver::Solve().

Vec LinearSystem::Solve ( Vec  lhsGuess = NULL  ) 

Solve the linear system.

Parameters:
lhsGuess an optional initial guess for the solution (defaults to NULL)

Todo:
never tested in linalg component

Todo:
Should it be compulsory for the caller to supply this and manage the memory?

Definition at line 539 of file LinearSystem.cpp.

References GenericEventHandler< 11, HeartEventHandler >::BeginEvent(), GenericEventHandler< 11, HeartEventHandler >::EndEvent(), mKspIsSetup, mKspSolver, mKspType, mLhsMatrix, mMatNullSpace, mMatrixIsConstant, mNonZerosUsed, mpBlockDiagonalPC, mPcType, mRhsVector, mSize, mTolerance, and mUseAbsoluteTolerance.

Referenced by SimpleNewtonNonlinearSolver::Solve(), and AbstractLinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::StaticSolve().

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

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

Parameters:
row the row index
value the value to set this entry

Definition at line 264 of file LinearSystem.cpp.

References mOwnershipRangeHi, mOwnershipRangeLo, and mRhsVector.

Referenced by AbstractNonlinearElasticityAssembler< DIM >::ApplyBoundaryConditions(), BoundaryConditionsContainer< ELEM_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletToLinearProblem(), and BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::FinaliseAssembleSystem().

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

References mOwnershipRangeHi, mOwnershipRangeLo, and mRhsVector.

unsigned LinearSystem::GetSize (  )  const

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

Todo:
Document this method and its parameters!
Parameters:
nullbasis 
numberOfBases 

Definition at line 390 of file LinearSystem.cpp.

References mMatNullSpace.

Referenced by BidomainDg0Assembler< ELEMENT_DIM, SPACE_DIM >::FinaliseAssembleSystem().

Vec & LinearSystem::rGetRhsVector (  ) 

Vec LinearSystem::GetRhsVector (  )  const

Get access to the RHS vector for archiving

Definition at line 445 of file LinearSystem.cpp.

References mRhsVector.

Mat & LinearSystem::rGetLhsMatrix (  ) 

Mat LinearSystem::GetLhsMatrix (  )  const

Get access to the LHS matrix for archiving

Definition at line 455 of file LinearSystem.cpp.

References mLhsMatrix.

Vec & LinearSystem::rGetDirichletBoundaryConditionsVector (  ) 

Get access to the dirichlet boundary conditions vector.

Should only be used by the BoundaryConditionsContainer.

Definition at line 460 of file LinearSystem.cpp.

References mDirichletBoundaryConditionsVector.

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

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

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

References mOwnershipRangeHi, and mOwnershipRangeLo.

double LinearSystem::GetMatrixElement ( PetscInt  row,
PetscInt  col 
)

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

References mLhsMatrix, mOwnershipRangeHi, and mOwnershipRangeLo.

Referenced by ZeroMatrixColumn().

double LinearSystem::GetRhsVectorElement ( PetscInt  row  ) 

Return an element of the RHS vector. May only be called for elements you own.

Parameters:
row the row index

Definition at line 416 of file LinearSystem.cpp.

References mOwnershipRangeHi, mOwnershipRangeLo, and mRhsVector.

unsigned LinearSystem::GetNumIterations (  )  const

Return the number of iterations taken by the last Solve()

Definition at line 429 of file LinearSystem.cpp.

References mKspIsSetup, and mKspSolver.

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 438 of file LinearSystem.hpp.

References mLhsMatrix, mOwnershipRangeHi, and mOwnershipRangeLo.

Referenced by NonlinearElasticityAssembler< DIM >::AssembleSystem(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleSystem().

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 501 of file LinearSystem.hpp.

References mOwnershipRangeHi, mOwnershipRangeLo, and mRhsVector.

Referenced by NonlinearElasticityAssembler< DIM >::AssembleSystem(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, NON_HEART, CONCRETE >::AssembleSystem().


Friends And Related Function Documentation

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

Needed for serialization.

Definition at line 104 of file LinearSystem.hpp.


Member Data Documentation

Mat LinearSystem::mLhsMatrix [private]

Vec LinearSystem::mRhsVector [private]

PetscInt LinearSystem::mSize [private]

The size of the linear system.

Definition at line 66 of file LinearSystem.hpp.

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

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 72 of file LinearSystem.hpp.

Referenced by AddLhsMultipleValues(), AddRhsMultipleValues(), AddToMatrixElement(), AddToRhsVectorElement(), GetMatrixElement(), GetOwnershipRange(), GetRhsVectorElement(), LinearSystem(), SetMatrixElement(), SetMatrixRow(), SetRhsVectorElement(), SetupVectorAndMatrix(), ZeroMatrixColumn(), and ZeroRhsVector().

PetscInt LinearSystem::mOwnershipRangeHi [private]

MatNullSpace LinearSystem::mMatNullSpace [private]

PETSc null matrix.

Definition at line 75 of file LinearSystem.hpp.

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

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

Definition at line 78 of file LinearSystem.hpp.

Referenced by ~LinearSystem().

KSP LinearSystem::mKspSolver [private]

The PETSc linear solver object

Definition at line 80 of file LinearSystem.hpp.

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

bool LinearSystem::mKspIsSetup [private]

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

Definition at line 81 of file LinearSystem.hpp.

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

double LinearSystem::mNonZerosUsed [private]

Yes, it really is stored as a double.

Definition at line 82 of file LinearSystem.hpp.

Referenced by serialize(), and Solve().

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

Definition at line 83 of file LinearSystem.hpp.

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

double LinearSystem::mTolerance [private]

absolute or relative tolerance of the KSP solver

Definition at line 84 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 89 of file LinearSystem.hpp.

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

std::string LinearSystem::mKspType [private]

KSP solver type (see PETSc KSPSetType() )

Definition at line 90 of file LinearSystem.hpp.

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

std::string LinearSystem::mPcType [private]

Preconditioner type (see PETSc PCSetType() )

Definition at line 91 of file LinearSystem.hpp.

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

Storage for efficient application of Dirichlet BCs, see AbstractBoundaryConditionsContainer

Definition at line 93 of file LinearSystem.hpp.

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

Stores a pointer to a purpose-build preconditioner

Definition at line 96 of file LinearSystem.hpp.

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


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

Generated on Tue Aug 4 16:11:21 2009 for Chaste by  doxygen 1.5.5