6 #ifndef CT_SQUAREMATRIX_H 7 #define CT_SQUAREMATRIX_H 38 int solve(doublereal* b,
size_t nrhs=1,
size_t ldb=0);
40 void resize(
size_t n,
size_t m, doublereal v = 0.0);
45 virtual void mult(
const doublereal* b, doublereal* prod)
const;
47 virtual void leftMult(
const doublereal*
const b, doublereal*
const prod)
const;
53 virtual doublereal
rcond(doublereal a1norm);
55 virtual doublereal
oneNorm()
const;
84 virtual doublereal
operator()(
size_t i,
size_t j)
const {
88 virtual size_t nRows()
const;
104 virtual vector_fp::iterator
begin();
105 virtual vector_fp::const_iterator
begin()
const;
107 virtual doublereal*
const* colPts();
109 virtual size_t checkRows(doublereal& valueSmall)
const;
110 virtual size_t checkColumns(doublereal& valueSmall)
const;
virtual size_t checkRows(doublereal &valueSmall) const
Check to see if we have any zero rows in the Jacobian.
virtual doublereal rcond(doublereal a1norm)
Returns an estimate of the inverse of the condition number for the matrix.
vector_fp work
Work vector for QR algorithm.
virtual size_t checkColumns(doublereal &valueSmall) const
Check to see if we have any zero columns in the Jacobian.
virtual doublereal operator()(size_t i, size_t j) const
Constant Index into the (i,j) element.
std::vector< int > vector_int
Vector of ints.
size_t nRowsAndStruct(size_t *const iStruct=0) const
Return the size and structure of the matrix.
vector_fp tau
Work vector for QR algorithm.
virtual vector_fp::iterator begin()
Return an iterator pointing to the first element.
int solve(doublereal *b, size_t nrhs=1, size_t ldb=0)
Solves the Ax = b system returning x in the b spot.
virtual size_t nRows() const
Return the number of rows in the matrix.
virtual doublereal & operator()(size_t i, size_t j)
Index into the (i,j) element.
void zero()
Zero the matrix.
A class for full (non-sparse) matrices with Fortran-compatible data storage.
int factor()
Factors the A matrix, overwriting A.
virtual int factorAlgorithm() const
Returns the factor algorithm used.
virtual void useFactorAlgorithm(int fAlgorithm)
Change the way the matrix is factored.
int solveQR(doublereal *b)
Solves the linear problem Ax=b using the QR algorithm returning x in the b spot.
virtual GeneralMatrix * duplMyselfAsGeneralMatrix() const
Duplicator member function.
doublereal & operator()(size_t i, size_t j)
Allows setting elements using the syntax A(i,j) = x.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Headers for the DenseMatrix object, which deals with dense rectangular matrices and description of th...
vector_int iwork_
Integer work vector for QR algorithms.
virtual doublereal oneNorm() const
Calculate the one norm of the matrix.
virtual doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are assumed to be contiguous in memory.
Declarations for the class GeneralMatrix which is a virtual base class for matrices handled by solver...
virtual void leftMult(const doublereal *const b, doublereal *const prod) const
Left-multiply the matrix by transpose(b), and write the result to prod.
virtual doublereal rcondQR()
Returns an estimate of the inverse of the condition number for the matrix.
int useQR_
Use the QR algorithm to factor and invert the matrix.
void setFactorFlag()
set the factored flag
Namespace for the Cantera kernel.
doublereal a1norm_
1-norm of the matrix.
SquareMatrix()
Base Constructor.
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
virtual int factorQR()
Factors the A matrix using the QR algorithm, overwriting A.