Cantera  3.1.0a1
BandMatrix Class Reference

A class for banded matrices, involving matrix inversion processes. More...

#include <BandMatrix.h>

Inheritance diagram for BandMatrix:
[legend]

Detailed Description

A class for banded matrices, involving matrix inversion processes.

The class is based upon the LAPACK banded storage matrix format.

An important issue with this class is that it stores both the original data and the LU factorization of the data. This means that the banded matrix typically will take up twice the room that it is expected to take.

QR factorizations of banded matrices are not included in the original LAPACK work. Add-ons are available. However, they are not included here. Instead we just use the stock LU decompositions.

This class is a derived class of the base class GeneralMatrix. However, within the oneD directory, the class is used as is, without reference to the GeneralMatrix base type.

Definition at line 36 of file BandMatrix.h.

Public Member Functions

 BandMatrix ()
 Base Constructor. More...
 
 BandMatrix (size_t n, size_t kl, size_t ku, double v=0.0)
 Creates a banded matrix and sets all elements to zero. More...
 
 BandMatrix (const BandMatrix &y)
 
BandMatrixoperator= (const BandMatrix &y)
 
void resize (size_t n, size_t kl, size_t ku, double v=0.0)
 Resize the matrix problem. More...
 
void bfill (double v=0.0)
 Fill or zero the matrix. More...
 
double & operator() (size_t i, size_t j) override
 Index into the (i,j) element. More...
 
double operator() (size_t i, size_t j) const override
 Constant Index into the (i,j) element. More...
 
double & value (size_t i, size_t j)
 Return a changeable reference to element (i,j). More...
 
double value (size_t i, size_t j) const
 Return the value of element (i,j). More...
 
size_t index (size_t i, size_t j) const
 Returns the location in the internal 1D array corresponding to the (i,j) element in the banded array. More...
 
double _value (size_t i, size_t j) const
 Return the value of the (i,j) element for (i,j) within the bandwidth. More...
 
size_t nRows () const override
 Return the number of rows in the matrix. More...
 
size_t nColumns () const
 Number of columns. More...
 
size_t nSubDiagonals () const
 Number of subdiagonals. More...
 
size_t nSuperDiagonals () const
 Number of superdiagonals. More...
 
size_t ldim () const
 Return the number of rows of storage needed for the band storage. More...
 
void mult (const double *b, double *prod) const override
 Multiply A*b and write result to prod. More...
 
void leftMult (const double *const b, double *const prod) const override
 Multiply b*A and write result to prod. More...
 
int factor () override
 Perform an LU decomposition, the LAPACK routine DGBTRF is used. More...
 
int solve (const double *const b, double *const x)
 Solve the matrix problem Ax = b. More...
 
int solve (double *b, size_t nrhs=1, size_t ldb=0) override
 Solve the matrix problem Ax = b. More...
 
void zero () override
 Zero the matrix elements. More...
 
double rcond (double a1norm) override
 Returns an estimate of the inverse of the condition number for the matrix. More...
 
int factorAlgorithm () const override
 Returns the factor algorithm used. More...
 
double oneNorm () const override
 Returns the one norm of the matrix. More...
 
double * ptrColumn (size_t j) override
 Return a pointer to the top of column j. More...
 
double *const * colPts () override
 Return a vector of const pointers to the columns. More...
 
size_t checkRows (double &valueSmall) const override
 Check to see if we have any zero rows in the Jacobian. More...
 
size_t checkColumns (double &valueSmall) const override
 Check to see if we have any zero columns in the Jacobian. More...
 
int info () const
 LAPACK "info" flag after last factor/solve operation. More...
 
- Public Member Functions inherited from GeneralMatrix
 GeneralMatrix ()=default
 Base Constructor. More...
 
virtual int factorQR ()
 Factors the A matrix using the QR algorithm, overwriting A. More...
 
virtual double rcondQR ()
 Returns an estimate of the inverse of the condition number for the matrix. More...
 
virtual void useFactorAlgorithm (int fAlgorithm)
 Change the way the matrix is factored. More...
 
virtual void clearFactorFlag ()
 clear the factored flag More...
 
virtual bool factored () const
 true if the current factorization is up to date with the matrix More...
 

Protected Attributes

vector< double > data
 Matrix data. More...
 
vector< double > ludata
 Factorized data. More...
 
size_t m_n = 0
 Number of rows and columns of the matrix. More...
 
size_t m_kl = 0
 Number of subdiagonals of the matrix. More...
 
size_t m_ku = 0
 Number of super diagonals of the matrix. More...
 
double m_zero = 0
 value of zero More...
 
unique_ptr< PivData > m_ipiv
 Pivot vector. More...
 
vector< double * > m_colPtrs
 Vector of column pointers. More...
 
vector< double * > m_lu_col_ptrs
 
vector< int > iwork_
 Extra work array needed - size = n. More...
 
vector< double > work_
 Extra dp work array needed - size = 3n. More...
 
int m_info = 0
 
- Protected Attributes inherited from GeneralMatrix
int m_factored = false
 Indicates whether the matrix is factored. More...
 

Constructor & Destructor Documentation

◆ BandMatrix() [1/2]

Base Constructor.

Create an 0 by 0 matrix, and initialize all elements to 0.

Definition at line 34 of file BandMatrix.cpp.

◆ BandMatrix() [2/2]

BandMatrix ( size_t  n,
size_t  kl,
size_t  ku,
double  v = 0.0 
)

Creates a banded matrix and sets all elements to zero.

Create an n by n banded matrix, and initialize all elements to v.

Parameters
nsize of the square matrix
klband size on the lower portion of the matrix
kuband size on the upper portion of the matrix
vinitial value of all matrix components.

Definition at line 44 of file BandMatrix.cpp.

Member Function Documentation

◆ resize()

void resize ( size_t  n,
size_t  kl,
size_t  ku,
double  v = 0.0 
)

Resize the matrix problem.

All data is lost

Parameters
nsize of the square matrix
klband size on the lower portion of the matrix
kuband size on the upper portion of the matrix
vinitial value of all matrix components.

Definition at line 107 of file BandMatrix.cpp.

◆ bfill()

void bfill ( double  v = 0.0)

Fill or zero the matrix.

Parameters
vFill value, defaults to zero.

Definition at line 126 of file BandMatrix.cpp.

◆ operator()() [1/2]

double & operator() ( size_t  i,
size_t  j 
)
overridevirtual

Index into the (i,j) element.

Parameters
irow
jcolumn
Returns
a changeable reference to the matrix entry

Implements GeneralMatrix.

Definition at line 138 of file BandMatrix.cpp.

◆ operator()() [2/2]

double operator() ( size_t  i,
size_t  j 
) const
overridevirtual

Constant Index into the (i,j) element.

Parameters
irow
jcolumn
Returns
an unchangeable reference to the matrix entry

Implements GeneralMatrix.

Definition at line 143 of file BandMatrix.cpp.

◆ value() [1/2]

double & value ( size_t  i,
size_t  j 
)

Return a changeable reference to element (i,j).

Since this method may alter the element value, it may need to be refactored, so the flag m_factored is set to false.

Parameters
irow
jcolumn
Returns
a reference to the value of the matrix entry

Definition at line 148 of file BandMatrix.cpp.

◆ value() [2/2]

double value ( size_t  i,
size_t  j 
) const

Return the value of element (i,j).

This method does not alter the array.

Parameters
irow
jcolumn
Returns
the value of the matrix entry

Definition at line 157 of file BandMatrix.cpp.

◆ index()

size_t index ( size_t  i,
size_t  j 
) const

Returns the location in the internal 1D array corresponding to the (i,j) element in the banded array.

Parameters
irow
jcolumn
Returns
the index of the matrix entry

Definition at line 165 of file BandMatrix.cpp.

◆ _value()

double _value ( size_t  i,
size_t  j 
) const

Return the value of the (i,j) element for (i,j) within the bandwidth.

For efficiency, this method does not check that (i,j) are within the bandwidth; it is up to the calling program to insure that this is true.

Parameters
irow
jcolumn
Returns
the value of the matrix entry

Definition at line 170 of file BandMatrix.cpp.

◆ nRows()

size_t nRows ( ) const
overridevirtual

Return the number of rows in the matrix.

Implements GeneralMatrix.

Definition at line 175 of file BandMatrix.cpp.

◆ nColumns()

size_t nColumns ( ) const

Number of columns.

Definition at line 180 of file BandMatrix.cpp.

◆ nSubDiagonals()

size_t nSubDiagonals ( ) const

Number of subdiagonals.

Definition at line 185 of file BandMatrix.cpp.

◆ nSuperDiagonals()

size_t nSuperDiagonals ( ) const

Number of superdiagonals.

Definition at line 190 of file BandMatrix.cpp.

◆ ldim()

size_t ldim ( ) const

Return the number of rows of storage needed for the band storage.

Definition at line 195 of file BandMatrix.cpp.

◆ mult()

void mult ( const double *  b,
double *  prod 
) const
overridevirtual

Multiply A*b and write result to prod.

Implements GeneralMatrix.

Definition at line 200 of file BandMatrix.cpp.

◆ leftMult()

void leftMult ( const double *const  b,
double *const  prod 
) const
overridevirtual

Multiply b*A and write result to prod.

Parameters
bVector to do the lh multiplication
prodOUTPUT vector to receive the result

Implements GeneralMatrix.

Definition at line 213 of file BandMatrix.cpp.

◆ factor()

int factor ( )
overridevirtual

Perform an LU decomposition, the LAPACK routine DGBTRF is used.

The factorization is saved in ludata.

Returns
a success flag. 0 indicates a success; ~0 indicates some error occurred, see the LAPACK documentation

Implements GeneralMatrix.

Definition at line 226 of file BandMatrix.cpp.

◆ solve() [1/2]

int solve ( const double *const  b,
double *const  x 
)

Solve the matrix problem Ax = b.

Parameters
bINPUT RHS of the problem
xOUTPUT solution to the problem
Returns
a success flag. 0 indicates a success; ~0 indicates some error occurred, see the LAPACK documentation

Definition at line 253 of file BandMatrix.cpp.

◆ solve() [2/2]

int solve ( double *  b,
size_t  nrhs = 1,
size_t  ldb = 0 
)
overridevirtual

Solve the matrix problem Ax = b.

Parameters
bINPUT RHS of the problem OUTPUT solution to the problem
nrhsNumber of right hand sides to solve
ldbLeading dimension of b. Default is nColumns()
Returns
a success flag. 0 indicates a success; ~0 indicates some error occurred, see the LAPACK documentation

Implements GeneralMatrix.

Definition at line 259 of file BandMatrix.cpp.

◆ zero()

void zero ( )
overridevirtual

Zero the matrix elements.

Implements GeneralMatrix.

Definition at line 132 of file BandMatrix.cpp.

◆ rcond()

double rcond ( double  a1norm)
overridevirtual

Returns an estimate of the inverse of the condition number for the matrix.

The matrix must have been previously factored using the LU algorithm

Parameters
a1normNorm of the matrix
Returns
the inverse of the condition number

Implements GeneralMatrix.

Definition at line 305 of file BandMatrix.cpp.

◆ factorAlgorithm()

int factorAlgorithm ( ) const
overridevirtual

Returns the factor algorithm used.

This method will always return 0 (LU) for band matrices.

Implements GeneralMatrix.

Definition at line 328 of file BandMatrix.cpp.

◆ oneNorm()

double oneNorm ( ) const
overridevirtual

Returns the one norm of the matrix.

Implements GeneralMatrix.

Definition at line 333 of file BandMatrix.cpp.

◆ ptrColumn()

double * ptrColumn ( size_t  j)
overridevirtual

Return a pointer to the top of column j.

Column values are assumed to be contiguous in memory (LAPACK band matrix structure)

On entry, the matrix A in band storage, in rows KL+1 to 2*KL+KU+1; rows 1 to KL of the array need not be set. The j-th column of A is stored in the j-th column of the array AB as follows:

AB(KL + KU + 1 + i - j,j) = A(i,j) for max(1, j - KU) <= i <= min(m, j + KL)

This routine returns the position of AB(1,j) (fortran-1 indexing) in the above format

So to address the (i,j) position, you use the following indexing:

double *colP_j = matrix.ptrColumn(j);
double a_i_j = colP_j[kl + ku + i - j];
Parameters
jValue of the column
Returns
a pointer to the top of the column

Implements GeneralMatrix.

Definition at line 392 of file BandMatrix.cpp.

◆ colPts()

double *const * colPts ( )
overridevirtual

Return a vector of const pointers to the columns.

Note the value of the pointers are protected by their being const. However, the value of the matrix is open to being changed.

Returns
a vector of pointers to the top of the columns of the matrices.

Implements GeneralMatrix.

Definition at line 397 of file BandMatrix.cpp.

◆ checkRows()

size_t checkRows ( double &  valueSmall) const
overridevirtual

Check to see if we have any zero rows in the Jacobian.

This utility routine checks to see if any rows are zero. The smallest row is returned along with the largest coefficient in that row

Parameters
valueSmallOUTPUT value of the largest coefficient in the smallest row
Returns
index of the row that is most nearly zero

Implements GeneralMatrix.

Definition at line 348 of file BandMatrix.cpp.

◆ checkColumns()

size_t checkColumns ( double &  valueSmall) const
overridevirtual

Check to see if we have any zero columns in the Jacobian.

This utility routine checks to see if any columns are zero. The smallest column is returned along with the largest coefficient in that column

Parameters
valueSmallOUTPUT value of the largest coefficient in the smallest column
Returns
index of the column that is most nearly zero

Implements GeneralMatrix.

Definition at line 370 of file BandMatrix.cpp.

◆ info()

int info ( ) const
inline

LAPACK "info" flag after last factor/solve operation.

Definition at line 239 of file BandMatrix.h.

Member Data Documentation

◆ data

vector<double> data
protected

Matrix data.

Definition at line 243 of file BandMatrix.h.

◆ ludata

vector<double> ludata
protected

Factorized data.

Definition at line 246 of file BandMatrix.h.

◆ m_n

size_t m_n = 0
protected

Number of rows and columns of the matrix.

Definition at line 249 of file BandMatrix.h.

◆ m_kl

size_t m_kl = 0
protected

Number of subdiagonals of the matrix.

Definition at line 252 of file BandMatrix.h.

◆ m_ku

size_t m_ku = 0
protected

Number of super diagonals of the matrix.

Definition at line 255 of file BandMatrix.h.

◆ m_zero

double m_zero = 0
protected

value of zero

Definition at line 258 of file BandMatrix.h.

◆ m_ipiv

unique_ptr<PivData> m_ipiv
protected

Pivot vector.

Definition at line 263 of file BandMatrix.h.

◆ m_colPtrs

vector<double*> m_colPtrs
protected

Vector of column pointers.

Definition at line 266 of file BandMatrix.h.

◆ iwork_

vector<int> iwork_
protected

Extra work array needed - size = n.

Definition at line 270 of file BandMatrix.h.

◆ work_

vector<double> work_
protected

Extra dp work array needed - size = 3n.

Definition at line 273 of file BandMatrix.h.


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