Cantera  2.1.2
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);
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  virtual void clearFactorFlag();
73 
74  //! set the factored flag
75  void setFactorFlag();
76 
77  virtual bool factored() const;
78  virtual void useFactorAlgorithm(int fAlgorithm);
79 
80  //! Returns the factor algorithm used
81  /*!
82  * 0 LU decomposition
83  * 1 QR decomposition
84  *
85  * This routine will always return 0
86  */
87  virtual int factorAlgorithm() const;
88 
89  virtual doublereal* ptrColumn(size_t j);
90 
91  virtual doublereal& operator()(size_t i, size_t j) {
92  return Array2D::operator()(i, j);
93  }
94 
95  virtual void copyData(const GeneralMatrix& y);
96 
97  virtual doublereal operator()(size_t i, size_t j) const {
98  return Array2D::operator()(i, j);
99  }
100 
101  virtual size_t nRows() const;
102 
103  //! Return the size and structure of the matrix
104  /*!
105  * This is inherited from GeneralMatrix
106  *
107  * @param iStruct OUTPUT Pointer to a vector of ints that describe the structure of the matrix.
108  * not used
109  *
110  * @return returns the number of rows and columns in the matrix.
111  */
112  size_t nRowsAndStruct(size_t* const iStruct = 0) const;
113 
114  virtual GeneralMatrix* duplMyselfAsGeneralMatrix() const;
115 
116  virtual vector_fp::iterator begin();
117  virtual vector_fp::const_iterator begin() const;
118 
119  virtual doublereal* const* colPts();
120 
121  virtual size_t checkRows(doublereal& valueSmall) const;
122  virtual size_t checkColumns(doublereal& valueSmall) const;
123 
124 protected:
125  //! the factor flag
127 
128 public:
129  //! Work vector for QR algorithm
131 
132  //! Work vector for QR algorithm
134 
135  //! Integer work vector for QR algorithms
136  std::vector<int> iwork_;
137 protected:
138  //! 1-norm of the matrix. This is determined immediately before every factorization
139  doublereal a1norm_;
140 
141  //! Use the QR algorithm to factor and invert the matrix
142  int useQR_;
143 };
144 }
145 
146 #endif
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.
Definition: SquareMatrix.h:133
Generic matrix.
Definition: GeneralMatrix.h:23
SquareMatrix & operator=(const SquareMatrix &right)
Assignment operator.
int m_factored
the factor flag
Definition: SquareMatrix.h:126
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.
Definition: SquareMatrix.h:97
vector_fp tau
Work vector for QR algorithm.
Definition: SquareMatrix.h:130
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.
Definition: SquareMatrix.h:91
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:274
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:165
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.
Definition: SquareMatrix.h:136
int useQR_
Use the QR algorithm to factor and invert the matrix.
Definition: SquareMatrix.h:142
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:139
SquareMatrix()
Base Constructor.
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition: DenseMatrix.h:70
virtual int factorQR()
Factors the A matrix using the QR algorithm, overwriting A.
virtual void clearFactorFlag()
clear the factored flag