Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SquareMatrix.h
Go to the documentation of this file.
1 /**
2  * @file SquareMatrix.h
3  * Dense, Square (not sparse) matrices.
4  */
5 
6 /*
7  * Copyright 2004 Sandia Corporation. Under the terms of Contract
8  * DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government
9  * retains certain rights in this software.
10  * See file License.txt for licensing information.
11  */
12 
13 #ifndef CT_SQUAREMATRIX_H
14 #define CT_SQUAREMATRIX_H
15 
16 #include "DenseMatrix.h"
17 #include "GeneralMatrix.h"
18 
19 namespace Cantera
20 {
21 
22 /**
23  * A class for full (non-sparse) matrices with Fortran-compatible
24  * data storage. Adds matrix inversion operations to this class from DenseMatrix.
25  */
26 class SquareMatrix: public DenseMatrix, public GeneralMatrix
27 {
28 public:
29  //! Base Constructor.
30  SquareMatrix();
31 
32  //! Constructor.
33  /*!
34  * Create an \c n by \c n matrix, and initialize all elements to \c v.
35  *
36  * @param n size of the square matrix
37  * @param v initial value of all matrix components.
38  */
39  SquareMatrix(size_t n, doublereal v = 0.0);
40 
41  //! Copy Constructor
42  SquareMatrix(const SquareMatrix& right);
43 
44  //! Assignment operator
45  SquareMatrix& operator=(const SquareMatrix& right);
46 
47  int solve(doublereal* b, size_t nrhs=1, size_t ldb=0);
48 
49  void resize(size_t n, size_t m, doublereal v = 0.0);
50 
51  //! Zero the matrix
52  void zero();
53 
54  virtual void mult(const doublereal* b, doublereal* prod) const;
55  virtual void mult(const DenseMatrix& b, DenseMatrix& prod) const;
56  virtual void leftMult(const doublereal* const b, doublereal* const prod) const;
57 
58  int factor();
59  virtual int factorQR();
60 
61  virtual doublereal rcondQR();
62  virtual doublereal rcond(doublereal a1norm);
63 
64  virtual doublereal oneNorm() const;
65 
66  //! Solves the linear problem Ax=b using the QR algorithm returning x in the b spot
67  /*!
68  * @param b RHS to be solved.
69  */
70  int solveQR(doublereal* b);
71 
72  //! set the factored flag
73  void setFactorFlag();
74 
75  virtual void useFactorAlgorithm(int fAlgorithm);
76 
77  //! Returns the factor algorithm used
78  /*!
79  * 0 LU decomposition
80  * 1 QR decomposition
81  *
82  * This routine will always return 0
83  */
84  virtual int factorAlgorithm() const;
85 
86  virtual doublereal* ptrColumn(size_t j);
87 
88  virtual doublereal& operator()(size_t i, size_t j) {
89  return Array2D::operator()(i, j);
90  }
91 
92  //! @deprecated To be removed after Cantera 2.2.
93  virtual void copyData(const GeneralMatrix& y);
94 
95  virtual doublereal operator()(size_t i, size_t j) const {
96  return Array2D::operator()(i, j);
97  }
98 
99  virtual size_t nRows() const;
100 
101  //! Return the size and structure of the matrix
102  /*!
103  * This is inherited from GeneralMatrix
104  *
105  * @param iStruct OUTPUT Pointer to a vector of ints that describe the structure of the matrix.
106  * not used
107  *
108  * @return returns the number of rows and columns in the matrix.
109  */
110  size_t nRowsAndStruct(size_t* const iStruct = 0) const;
111 
112  virtual GeneralMatrix* duplMyselfAsGeneralMatrix() const;
113 
114  virtual vector_fp::iterator begin();
115  virtual vector_fp::const_iterator begin() const;
116 
117  virtual doublereal* const* colPts();
118 
119  virtual size_t checkRows(doublereal& valueSmall) const;
120  virtual size_t checkColumns(doublereal& valueSmall) const;
121 
122  //! Work vector for QR algorithm
124 
125  //! Work vector for QR algorithm
127 
128  //! Integer work vector for QR algorithms
129  std::vector<int> iwork_;
130 protected:
131  //! 1-norm of the matrix. This is determined immediately before every factorization
132  doublereal a1norm_;
133 
134  //! Use the QR algorithm to factor and invert the matrix
135  int useQR_;
136 };
137 }
138 
139 #endif
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.
Definition: SquareMatrix.h:126
Generic matrix.
Definition: GeneralMatrix.h:24
SquareMatrix & operator=(const SquareMatrix &right)
Assignment operator.
virtual GeneralMatrix * duplMyselfAsGeneralMatrix() const
Duplicator member function.
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.
Definition: SquareMatrix.h:95
vector_fp tau
Work vector for QR algorithm.
Definition: SquareMatrix.h:123
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 doublereal & operator()(size_t i, size_t j)
Index into the (i,j) element.
Definition: SquareMatrix.h:88
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.
Definition: SquareMatrix.h:26
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.
Definition: Array.h:269
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.
Definition: ct_defs.h:157
virtual void copyData(const GeneralMatrix &y)
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.
Definition: SquareMatrix.h:129
int useQR_
Use the QR algorithm to factor and invert the matrix.
Definition: SquareMatrix.h:135
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
Definition: SquareMatrix.h:132
SquareMatrix()
Base Constructor.
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition: DenseMatrix.h:71
virtual int factorQR()
Factors the A matrix using the QR algorithm, overwriting A.