Cantera
2.5.1
|
A class for banded matrices, involving matrix inversion processes. More...
#include <BandMatrix.h>
Public Member Functions | |
BandMatrix () | |
Base Constructor. More... | |
BandMatrix (size_t n, size_t kl, size_t ku, doublereal v=0.0) | |
Creates a banded matrix and sets all elements to zero. More... | |
BandMatrix (const BandMatrix &y) | |
BandMatrix & | operator= (const BandMatrix &y) |
void | resize (size_t n, size_t kl, size_t ku, doublereal v=0.0) |
Resize the matrix problem. More... | |
void | bfill (doublereal v=0.0) |
Fill or zero the matrix. More... | |
doublereal & | operator() (size_t i, size_t j) |
Index into the (i,j) element. More... | |
doublereal | operator() (size_t i, size_t j) const |
Constant Index into the (i,j) element. More... | |
doublereal & | value (size_t i, size_t j) |
Return a changeable reference to element (i,j). More... | |
doublereal | 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... | |
doublereal | _value (size_t i, size_t j) const |
Return the value of the (i,j) element for (i,j) within the bandwidth. More... | |
virtual size_t | nRows () const |
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... | |
virtual void | mult (const doublereal *b, doublereal *prod) const |
Multiply A*b and write result to prod . More... | |
virtual void | leftMult (const doublereal *const b, doublereal *const prod) const |
Multiply b*A and write result to prod. More... | |
int | factor () |
Perform an LU decomposition, the LAPACK routine DGBTRF is used. More... | |
int | solve (const doublereal *const b, doublereal *const x) |
Solve the matrix problem Ax = b. More... | |
int | solve (doublereal *b, size_t nrhs=1, size_t ldb=0) |
Solve the matrix problem Ax = b. More... | |
virtual vector_fp::iterator | begin () |
Returns an iterator for the start of the band storage data. More... | |
vector_fp::iterator | end () |
Returns an iterator for the end of the band storage data. More... | |
vector_fp::const_iterator | begin () const |
Returns a const iterator for the start of the band storage data. More... | |
vector_fp::const_iterator | end () const |
Returns a const iterator for the end of the band storage data. More... | |
virtual void | zero () |
Zero the matrix elements. More... | |
virtual doublereal | rcond (doublereal a1norm) |
Returns an estimate of the inverse of the condition number for the matrix. More... | |
virtual int | factorAlgorithm () const |
Returns the factor algorithm used. More... | |
virtual doublereal | oneNorm () const |
Returns the one norm of the matrix. More... | |
virtual doublereal * | ptrColumn (size_t j) |
Return a pointer to the top of column j. More... | |
virtual doublereal *const * | colPts () |
Return a vector of const pointers to the columns. More... | |
virtual size_t | checkRows (doublereal &valueSmall) const |
Check to see if we have any zero rows in the Jacobian. More... | |
virtual size_t | checkColumns (doublereal &valueSmall) const |
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 () | |
Base Constructor. More... | |
virtual int | factorQR () |
Factors the A matrix using the QR algorithm, overwriting A. More... | |
virtual doublereal | 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_fp | data |
Matrix data. More... | |
vector_fp | ludata |
Factorized data. More... | |
size_t | m_n |
Number of rows and columns of the matrix. More... | |
size_t | m_kl |
Number of subdiagonals of the matrix. More... | |
size_t | m_ku |
Number of super diagonals of the matrix. More... | |
doublereal | m_zero |
value of zero More... | |
std::unique_ptr< PivData > | m_ipiv |
Pivot vector. More... | |
std::vector< doublereal * > | m_colPtrs |
Vector of column pointers. More... | |
std::vector< double * > | m_lu_col_ptrs |
vector_int | iwork_ |
Extra work array needed - size = n. More... | |
vector_fp | work_ |
Extra dp work array needed - size = 3n. More... | |
int | m_info |
Protected Attributes inherited from GeneralMatrix | |
int | m_factored |
Indicates whether the matrix is factored. More... | |
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 34 of file BandMatrix.h.
BandMatrix | ( | ) |
Base Constructor.
Create an 0
by 0
matrix, and initialize all elements to 0
.
Definition at line 51 of file BandMatrix.cpp.
BandMatrix | ( | size_t | n, |
size_t | kl, | ||
size_t | ku, | ||
doublereal | 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 66 of file BandMatrix.cpp.
References BandMatrix::data, BandMatrix::ludata, BandMatrix::m_colPtrs, BandMatrix::m_ipiv, and BandMatrix::m_n.
void resize | ( | size_t | n, |
size_t | kl, | ||
size_t | ku, | ||
doublereal | 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 135 of file BandMatrix.cpp.
References BandMatrix::data, BandMatrix::ludata, BandMatrix::m_colPtrs, GeneralMatrix::m_factored, BandMatrix::m_ipiv, BandMatrix::m_kl, BandMatrix::m_ku, and BandMatrix::m_n.
void bfill | ( | doublereal | v = 0.0 | ) |
Fill or zero the matrix.
v | Fill value, defaults to zero. |
Definition at line 154 of file BandMatrix.cpp.
References BandMatrix::data, and GeneralMatrix::m_factored.
|
virtual |
Index into the (i,j) element.
i | row |
j | column |
Implements GeneralMatrix.
Definition at line 166 of file BandMatrix.cpp.
References BandMatrix::value().
|
virtual |
Constant Index into the (i,j) element.
i | row |
j | column |
Implements GeneralMatrix.
Definition at line 171 of file BandMatrix.cpp.
References BandMatrix::value().
doublereal & 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 176 of file BandMatrix.cpp.
References BandMatrix::data, BandMatrix::index(), GeneralMatrix::m_factored, BandMatrix::m_kl, and BandMatrix::m_zero.
Referenced by BandMatrix::oneNorm(), and BandMatrix::operator()().
doublereal 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 185 of file BandMatrix.cpp.
References BandMatrix::data, BandMatrix::index(), and BandMatrix::m_kl.
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 193 of file BandMatrix.cpp.
References BandMatrix::m_kl, and BandMatrix::m_ku.
Referenced by BandMatrix::_value(), and BandMatrix::value().
doublereal _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 198 of file BandMatrix.cpp.
References BandMatrix::data, and BandMatrix::index().
Referenced by BandMatrix::checkColumns(), BandMatrix::checkRows(), BandMatrix::leftMult(), BandMatrix::mult(), and BandMatrix::oneNorm().
|
virtual |
Return the number of rows in the matrix.
Implements GeneralMatrix.
Definition at line 203 of file BandMatrix.cpp.
References BandMatrix::m_n.
Referenced by BandMatrix::factor(), and Cantera::operator<<().
size_t nColumns | ( | ) | const |
Number of columns.
Definition at line 208 of file BandMatrix.cpp.
References BandMatrix::m_n.
Referenced by BandMatrix::factor(), and BandMatrix::solve().
size_t nSubDiagonals | ( | ) | const |
Number of subdiagonals.
Definition at line 213 of file BandMatrix.cpp.
References BandMatrix::m_kl.
Referenced by BandMatrix::factor().
size_t nSuperDiagonals | ( | ) | const |
Number of superdiagonals.
Definition at line 218 of file BandMatrix.cpp.
References BandMatrix::m_ku.
Referenced by BandMatrix::factor().
size_t ldim | ( | ) | const |
Return the number of rows of storage needed for the band storage.
Definition at line 223 of file BandMatrix.cpp.
References BandMatrix::m_kl, and BandMatrix::m_ku.
Referenced by BandMatrix::factor().
|
virtual |
Multiply A*b and write result to prod
.
Implements GeneralMatrix.
Definition at line 228 of file BandMatrix.cpp.
References BandMatrix::_value(), BandMatrix::m_kl, BandMatrix::m_ku, and BandMatrix::m_n.
|
virtual |
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 241 of file BandMatrix.cpp.
References BandMatrix::_value(), BandMatrix::m_kl, BandMatrix::m_ku, and BandMatrix::m_n.
|
virtual |
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
The factorization is saved in ludata.
Implements GeneralMatrix.
Definition at line 254 of file BandMatrix.cpp.
References BandMatrix::data, BandMatrix::ldim(), BandMatrix::ludata, GeneralMatrix::m_factored, BandMatrix::m_ipiv, BandMatrix::nColumns(), BandMatrix::nRows(), BandMatrix::nSubDiagonals(), and BandMatrix::nSuperDiagonals().
Referenced by BandMatrix::solve().
int solve | ( | const doublereal *const | b, |
doublereal *const | x | ||
) |
Solve the matrix problem Ax = b.
b | INPUT RHS of the problem |
x | OUTPUT solution to the problem |
Definition at line 275 of file BandMatrix.cpp.
References BandMatrix::m_n.
Referenced by MultiNewton::step().
|
virtual |
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 281 of file BandMatrix.cpp.
References BandMatrix::factor(), GeneralMatrix::m_factored, and BandMatrix::nColumns().
|
virtual |
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 309 of file BandMatrix.cpp.
References BandMatrix::data, and GeneralMatrix::m_factored.
vector_fp::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 315 of file BandMatrix.cpp.
References BandMatrix::data, and GeneralMatrix::m_factored.
|
virtual |
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 321 of file BandMatrix.cpp.
References BandMatrix::data.
vector_fp::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 326 of file BandMatrix.cpp.
References BandMatrix::data.
|
virtual |
Zero the matrix elements.
Implements GeneralMatrix.
Definition at line 160 of file BandMatrix.cpp.
References BandMatrix::data, and GeneralMatrix::m_factored.
|
virtual |
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 343 of file BandMatrix.cpp.
References BandMatrix::iwork_, BandMatrix::ludata, GeneralMatrix::m_factored, BandMatrix::m_ipiv, BandMatrix::m_kl, BandMatrix::m_ku, BandMatrix::m_n, and BandMatrix::work_.
|
virtual |
Returns the factor algorithm used.
This method will always return 0 (LU) for band matrices.
Implements GeneralMatrix.
Definition at line 366 of file BandMatrix.cpp.
|
virtual |
Returns the one norm of the matrix.
Implements GeneralMatrix.
Definition at line 371 of file BandMatrix.cpp.
References BandMatrix::_value(), BandMatrix::m_kl, BandMatrix::m_ku, BandMatrix::m_n, and BandMatrix::value().
|
virtual |
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 430 of file BandMatrix.cpp.
References BandMatrix::m_colPtrs.
|
virtual |
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 435 of file BandMatrix.cpp.
References BandMatrix::m_colPtrs.
|
virtual |
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 386 of file BandMatrix.cpp.
References BandMatrix::_value(), BandMatrix::m_kl, BandMatrix::m_ku, BandMatrix::m_n, and Cantera::npos.
|
virtual |
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 408 of file BandMatrix.cpp.
References BandMatrix::_value(), BandMatrix::m_kl, BandMatrix::m_ku, BandMatrix::m_n, and Cantera::npos.
|
inline |
LAPACK "info" flag after last factor/solve operation.
Definition at line 261 of file BandMatrix.h.
Referenced by MultiNewton::step().
|
protected |
Matrix data.
Definition at line 265 of file BandMatrix.h.
Referenced by BandMatrix::_value(), BandMatrix::BandMatrix(), BandMatrix::begin(), BandMatrix::bfill(), BandMatrix::end(), BandMatrix::factor(), BandMatrix::resize(), BandMatrix::value(), and BandMatrix::zero().
|
protected |
Factorized data.
Definition at line 268 of file BandMatrix.h.
Referenced by BandMatrix::BandMatrix(), BandMatrix::factor(), BandMatrix::rcond(), and BandMatrix::resize().
|
protected |
Number of rows and columns of the matrix.
Definition at line 271 of file BandMatrix.h.
Referenced by BandMatrix::BandMatrix(), BandMatrix::checkColumns(), BandMatrix::checkRows(), BandMatrix::leftMult(), BandMatrix::mult(), BandMatrix::nColumns(), BandMatrix::nRows(), BandMatrix::oneNorm(), BandMatrix::rcond(), BandMatrix::resize(), and BandMatrix::solve().
|
protected |
Number of subdiagonals of the matrix.
Definition at line 274 of file BandMatrix.h.
Referenced by BandMatrix::checkColumns(), BandMatrix::checkRows(), BandMatrix::index(), BandMatrix::ldim(), BandMatrix::leftMult(), BandMatrix::mult(), BandMatrix::nSubDiagonals(), BandMatrix::oneNorm(), BandMatrix::rcond(), BandMatrix::resize(), and BandMatrix::value().
|
protected |
Number of super diagonals of the matrix.
Definition at line 277 of file BandMatrix.h.
Referenced by BandMatrix::checkColumns(), BandMatrix::checkRows(), BandMatrix::index(), BandMatrix::ldim(), BandMatrix::leftMult(), BandMatrix::mult(), BandMatrix::nSuperDiagonals(), BandMatrix::oneNorm(), BandMatrix::rcond(), and BandMatrix::resize().
|
protected |
|
protected |
Pivot vector.
Definition at line 285 of file BandMatrix.h.
Referenced by BandMatrix::BandMatrix(), BandMatrix::factor(), BandMatrix::rcond(), and BandMatrix::resize().
|
protected |
Vector of column pointers.
Definition at line 288 of file BandMatrix.h.
Referenced by BandMatrix::BandMatrix(), BandMatrix::colPts(), BandMatrix::ptrColumn(), and BandMatrix::resize().
|
protected |
Extra work array needed - size = n.
Definition at line 292 of file BandMatrix.h.
Referenced by BandMatrix::rcond().
|
protected |
Extra dp work array needed - size = 3n.
Definition at line 295 of file BandMatrix.h.
Referenced by BandMatrix::rcond().