13 #ifndef CT_SQUAREMATRIX_H
14 #define CT_SQUAREMATRIX_H
47 int solve(doublereal* b);
49 void resize(
size_t n,
size_t m, doublereal v = 0.0);
54 virtual void mult(
const doublereal* b, doublereal* prod)
const;
56 virtual void leftMult(
const doublereal*
const b, doublereal*
const prod)
const;
62 virtual doublereal
rcond(doublereal a1norm);
64 virtual doublereal
oneNorm()
const;
97 virtual doublereal
operator()(
size_t i,
size_t j)
const {
101 virtual size_t nRows()
const;
116 virtual vector_fp::iterator
begin();
117 virtual vector_fp::const_iterator
begin()
const;
119 virtual doublereal*
const* colPts();
121 virtual size_t checkRows(doublereal& valueSmall)
const;
122 virtual size_t checkColumns(doublereal& valueSmall)
const;
virtual bool factored() const
true if the current factorization is up to date with the matrix
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.
SquareMatrix & operator=(const SquareMatrix &right)
Assignment operator.
int m_factored
the factor flag
virtual GeneralMatrix * duplMyselfAsGeneralMatrix() const
Duplicator member function.
int solve(doublereal *b)
Solves the Ax = b system returning x in the b spot.
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.
vector_fp tau
Work vector for QR algorithm.
virtual vector_fp::iterator begin()
Return an iterator pointing to the first element.
virtual doublereal & operator()(size_t i, size_t j)
Index into the (i,j) element.
virtual void leftMult(const doublereal *const b, doublereal *const prod) const
Left-multiply the matrix by transpose(b), and write the result to prod.
void zero()
Zero the matrix.
virtual doublereal oneNorm() const
Calculate the one norm of the matrix.
A class for full (non-sparse) matrices with Fortran-compatible data storage.
int factor()
Factors the A matrix, overwriting A.
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.
doublereal & operator()(size_t i, size_t j)
Allows setting elements using the syntax A(i,j) = x.
virtual size_t checkRows(doublereal &valueSmall) const
Check to see if we have any zero rows in the jacobian.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
virtual void copyData(const GeneralMatrix &y)
Copy the data from one array into another without doing any checking.
Headers for the DenseMatrix object, which deals with dense rectangular matrices and description of th...
virtual size_t nRows() const
Return the number of rows in the matrix.
virtual int factorAlgorithm() const
Returns the factor algorithm used.
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 doublereal rcondQR()
Returns an estimate of the inverse of the condition number for the matrix.
std::vector< int > iwork_
Integer work vector for QR algorithms.
int useQR_
Use the QR algorithm to factor and invert the matrix.
size_t nRowsAndStruct(size_t *const iStruct=0) const
Return the size and structure of the matrix.
void setFactorFlag()
set the factored flag
doublereal a1norm_
1-norm of the matrix. This is determined immediately before every factorization
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.
virtual void clearFactorFlag()
clear the factored flag