Cantera
3.0.0
|
A class for banded matrices, involving matrix inversion processes. More...
#include <BandMatrix.h>
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. | |
BandMatrix (size_t n, size_t kl, size_t ku, double v=0.0) | |
Creates a banded matrix and sets all elements to zero. | |
BandMatrix (const BandMatrix &y) | |
BandMatrix & | operator= (const BandMatrix &y) |
void | resize (size_t n, size_t kl, size_t ku, double v=0.0) |
Resize the matrix problem. | |
void | bfill (double v=0.0) |
Fill or zero the matrix. | |
double & | operator() (size_t i, size_t j) override |
Index into the (i,j) element. | |
double | operator() (size_t i, size_t j) const override |
Constant Index into the (i,j) element. | |
double & | value (size_t i, size_t j) |
Return a changeable reference to element (i,j). | |
double | value (size_t i, size_t j) const |
Return the value of element (i,j). | |
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. | |
double | _value (size_t i, size_t j) const |
Return the value of the (i,j) element for (i,j) within the bandwidth. | |
size_t | nRows () const override |
Return the number of rows in the matrix. | |
size_t | nColumns () const |
Number of columns. | |
size_t | nSubDiagonals () const |
Number of subdiagonals. | |
size_t | nSuperDiagonals () const |
Number of superdiagonals. | |
size_t | ldim () const |
Return the number of rows of storage needed for the band storage. | |
void | mult (const double *b, double *prod) const override |
Multiply A*b and write result to prod . | |
void | leftMult (const double *const b, double *const prod) const override |
Multiply b*A and write result to prod. | |
int | factor () override |
Perform an LU decomposition, the LAPACK routine DGBTRF is used. | |
int | solve (const double *const b, double *const x) |
Solve the matrix problem Ax = b. | |
int | solve (double *b, size_t nrhs=1, size_t ldb=0) override |
Solve the matrix problem Ax = b. | |
vector< double >::iterator | begin () override |
Returns an iterator for the start of the band storage data. | |
vector< double >::iterator | end () |
Returns an iterator for the end of the band storage data. | |
vector< double >::const_iterator | begin () const override |
Returns a const iterator for the start of the band storage data. | |
vector< double >::const_iterator | end () const |
Returns a const iterator for the end of the band storage data. | |
void | zero () override |
Zero the matrix elements. | |
double | rcond (double a1norm) override |
Returns an estimate of the inverse of the condition number for the matrix. | |
int | factorAlgorithm () const override |
Returns the factor algorithm used. | |
double | oneNorm () const override |
Returns the one norm of the matrix. | |
double * | ptrColumn (size_t j) override |
Return a pointer to the top of column j. | |
double *const * | colPts () override |
Return a vector of const pointers to the columns. | |
size_t | checkRows (double &valueSmall) const override |
Check to see if we have any zero rows in the Jacobian. | |
size_t | checkColumns (double &valueSmall) const override |
Check to see if we have any zero columns in the Jacobian. | |
int | info () const |
LAPACK "info" flag after last factor/solve operation. | |
Public Member Functions inherited from GeneralMatrix | |
GeneralMatrix ()=default | |
Base Constructor. | |
virtual void | zero ()=0 |
Zero the matrix elements. | |
virtual void | mult (const double *b, double *prod) const =0 |
Multiply A*b and write result to prod. | |
virtual void | leftMult (const double *const b, double *const prod) const =0 |
Multiply b*A and write result to prod. | |
virtual int | factor ()=0 |
Factors the A matrix, overwriting A. | |
virtual int | factorQR () |
Factors the A matrix using the QR algorithm, overwriting A. | |
virtual double | rcondQR () |
Returns an estimate of the inverse of the condition number for the matrix. | |
virtual double | rcond (double a1norm)=0 |
Returns an estimate of the inverse of the condition number for the matrix. | |
virtual void | useFactorAlgorithm (int fAlgorithm) |
Change the way the matrix is factored. | |
virtual int | factorAlgorithm () const =0 |
Return the factor algorithm used. | |
virtual double | oneNorm () const =0 |
Calculate the one norm of the matrix. | |
virtual size_t | nRows () const =0 |
Return the number of rows in the matrix. | |
virtual void | clearFactorFlag () |
clear the factored flag | |
virtual int | solve (double *b, size_t nrhs=1, size_t ldb=0)=0 |
Solves the Ax = b system returning x in the b spot. | |
virtual bool | factored () const |
true if the current factorization is up to date with the matrix | |
virtual double * | ptrColumn (size_t j)=0 |
Return a pointer to the top of column j, columns are assumed to be contiguous in memory. | |
virtual double & | operator() (size_t i, size_t j)=0 |
Index into the (i,j) element. | |
virtual double | operator() (size_t i, size_t j) const =0 |
Constant Index into the (i,j) element. | |
virtual vector< double >::iterator | begin ()=0 |
Return an iterator pointing to the first element. | |
virtual vector< double >::const_iterator | begin () const =0 |
Return a const iterator pointing to the first element. | |
virtual double *const * | colPts ()=0 |
Return a vector of const pointers to the columns. | |
virtual size_t | checkRows (double &valueSmall) const =0 |
Check to see if we have any zero rows in the Jacobian. | |
virtual size_t | checkColumns (double &valueSmall) const =0 |
Check to see if we have any zero columns in the Jacobian. | |
Protected Attributes | |
vector< double > | data |
Matrix data. | |
vector< double > | ludata |
Factorized data. | |
size_t | m_n = 0 |
Number of rows and columns of the matrix. | |
size_t | m_kl = 0 |
Number of subdiagonals of the matrix. | |
size_t | m_ku = 0 |
Number of super diagonals of the matrix. | |
double | m_zero = 0 |
value of zero | |
unique_ptr< PivData > | m_ipiv |
Pivot vector. | |
vector< double * > | m_colPtrs |
Vector of column pointers. | |
vector< double * > | m_lu_col_ptrs |
vector< int > | iwork_ |
Extra work array needed - size = n. | |
vector< double > | work_ |
Extra dp work array needed - size = 3n. | |
int | m_info = 0 |
Protected Attributes inherited from GeneralMatrix | |
int | m_factored = false |
Indicates whether the matrix is factored. | |
BandMatrix | ( | ) |
Base Constructor.
Create an 0
by 0
matrix, and initialize all elements to 0
.
Definition at line 34 of file BandMatrix.cpp.
~BandMatrix | ( | ) |
Definition at line 39 of file BandMatrix.cpp.
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
.
n | size of the square matrix |
kl | band size on the lower portion of the matrix |
ku | band size on the upper portion of the matrix |
v | initial value of all matrix components. |
Definition at line 44 of file BandMatrix.cpp.
BandMatrix | ( | const BandMatrix & | y | ) |
Definition at line 64 of file BandMatrix.cpp.
BandMatrix & operator= | ( | const BandMatrix & | y | ) |
Definition at line 84 of file BandMatrix.cpp.
void resize | ( | size_t | n, |
size_t | kl, | ||
size_t | ku, | ||
double | v = 0.0 |
||
) |
Resize the matrix problem.
All data is lost
n | size of the square matrix |
kl | band size on the lower portion of the matrix |
ku | band size on the upper portion of the matrix |
v | initial value of all matrix components. |
Definition at line 107 of file BandMatrix.cpp.
void bfill | ( | double | v = 0.0 | ) |
Fill or zero the matrix.
v | Fill value, defaults to zero. |
Definition at line 126 of file BandMatrix.cpp.
|
overridevirtual |
Index into the (i,j) element.
i | row |
j | column |
Implements GeneralMatrix.
Definition at line 138 of file BandMatrix.cpp.
|
overridevirtual |
Constant Index into the (i,j) element.
i | row |
j | column |
Implements GeneralMatrix.
Definition at line 143 of file BandMatrix.cpp.
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.
i | row |
j | column |
Definition at line 148 of file BandMatrix.cpp.
double value | ( | size_t | i, |
size_t | j | ||
) | const |
Return the value of element (i,j).
This method does not alter the array.
i | row |
j | column |
Definition at line 157 of file BandMatrix.cpp.
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.
i | row |
j | column |
Definition at line 165 of file BandMatrix.cpp.
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.
i | row |
j | column |
Definition at line 170 of file BandMatrix.cpp.
|
overridevirtual |
Return the number of rows in the matrix.
Implements GeneralMatrix.
Definition at line 175 of file BandMatrix.cpp.
size_t nColumns | ( | ) | const |
Number of columns.
Definition at line 180 of file BandMatrix.cpp.
size_t nSubDiagonals | ( | ) | const |
Number of subdiagonals.
Definition at line 185 of file BandMatrix.cpp.
size_t nSuperDiagonals | ( | ) | const |
Number of superdiagonals.
Definition at line 190 of file BandMatrix.cpp.
size_t ldim | ( | ) | const |
Return the number of rows of storage needed for the band storage.
Definition at line 195 of file BandMatrix.cpp.
|
overridevirtual |
Multiply A*b and write result to prod
.
Implements GeneralMatrix.
Definition at line 200 of file BandMatrix.cpp.
|
overridevirtual |
Multiply b*A and write result to prod.
b | Vector to do the lh multiplication |
prod | OUTPUT vector to receive the result |
Implements GeneralMatrix.
Definition at line 213 of file BandMatrix.cpp.
|
overridevirtual |
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
The factorization is saved in ludata.
Implements GeneralMatrix.
Definition at line 226 of file BandMatrix.cpp.
int solve | ( | const double *const | b, |
double *const | x | ||
) |
Solve the matrix problem Ax = b.
b | INPUT RHS of the problem |
x | OUTPUT solution to the problem |
Definition at line 253 of file BandMatrix.cpp.
|
overridevirtual |
Solve the matrix problem Ax = b.
b | INPUT RHS of the problem OUTPUT solution to the problem |
nrhs | Number of right hand sides to solve |
ldb | Leading dimension of b . Default is nColumns() |
Implements GeneralMatrix.
Definition at line 259 of file BandMatrix.cpp.
|
overridevirtual |
Returns an iterator for the start of the band storage data.
Iterator points to the beginning of the data, and it is changeable.
Implements GeneralMatrix.
Definition at line 293 of file BandMatrix.cpp.
vector< double >::iterator end | ( | ) |
Returns an iterator for the end of the band storage data.
Iterator points to the end of the data, and it is changeable.
Definition at line 300 of file BandMatrix.cpp.
|
overridevirtual |
Returns a const iterator for the start of the band storage data.
Iterator points to the beginning of the data, and it is not changeable.
Implements GeneralMatrix.
Definition at line 307 of file BandMatrix.cpp.
vector< double >::const_iterator end | ( | ) | const |
Returns a const iterator for the end of the band storage data.
Iterator points to the end of the data, and it is not changeable.
Definition at line 313 of file BandMatrix.cpp.
|
overridevirtual |
|
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
a1norm | Norm of the matrix |
Implements GeneralMatrix.
Definition at line 331 of file BandMatrix.cpp.
|
overridevirtual |
Returns the factor algorithm used.
This method will always return 0 (LU) for band matrices.
Implements GeneralMatrix.
Definition at line 354 of file BandMatrix.cpp.
|
overridevirtual |
Returns the one norm of the matrix.
Implements GeneralMatrix.
Definition at line 359 of file BandMatrix.cpp.
|
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];
j | Value of the column |
Implements GeneralMatrix.
Definition at line 418 of file BandMatrix.cpp.
|
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.
Implements GeneralMatrix.
Definition at line 423 of file BandMatrix.cpp.
|
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
valueSmall | OUTPUT value of the largest coefficient in the smallest row |
Implements GeneralMatrix.
Definition at line 374 of file BandMatrix.cpp.
|
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
valueSmall | OUTPUT value of the largest coefficient in the smallest column |
Implements GeneralMatrix.
Definition at line 396 of file BandMatrix.cpp.
|
inline |
LAPACK "info" flag after last factor/solve operation.
Definition at line 267 of file BandMatrix.h.
|
protected |
Matrix data.
Definition at line 271 of file BandMatrix.h.
|
protected |
Factorized data.
Definition at line 274 of file BandMatrix.h.
|
protected |
Number of rows and columns of the matrix.
Definition at line 277 of file BandMatrix.h.
|
protected |
Number of subdiagonals of the matrix.
Definition at line 280 of file BandMatrix.h.
|
protected |
Number of super diagonals of the matrix.
Definition at line 283 of file BandMatrix.h.
|
protected |
value of zero
Definition at line 286 of file BandMatrix.h.
|
protected |
Pivot vector.
Definition at line 291 of file BandMatrix.h.
|
protected |
Vector of column pointers.
Definition at line 294 of file BandMatrix.h.
|
protected |
Definition at line 295 of file BandMatrix.h.
|
protected |
Extra work array needed - size = n.
Definition at line 298 of file BandMatrix.h.
|
protected |
Extra dp work array needed - size = 3n.
Definition at line 301 of file BandMatrix.h.
|
protected |
Definition at line 303 of file BandMatrix.h.