Cantera  2.1.2
GeneralMatrix.h
Go to the documentation of this file.
1 /**
2  * @file GeneralMatrix.h
3  * Declarations for the class GeneralMatrix which is a virtual base class for matrices handled by solvers
4  * (see class \ref numerics and \link Cantera::GeneralMatrix GeneralMatrix\endlink).
5  */
6 
7 /*
8  * Copyright 2004 Sandia Corporation. Under the terms of Contract
9  * DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government
10  * retains certain rights in this software.
11  * See file License.txt for licensing information.
12  */
13 
14 #ifndef CT_GENERALMATRIX_H
15 #define CT_GENERALMATRIX_H
16 
17 #include "cantera/base/ct_defs.h"
18 
19 namespace Cantera
20 {
21 
22 //! Generic matrix
24 {
25 public:
26  //! Base Constructor
27  /*!
28  * @param matType Matrix type
29  * 0 full
30  * 1 banded
31  */
32  GeneralMatrix(int matType);
33 
34  //! Copy Constructor
35  GeneralMatrix(const GeneralMatrix& right);
36 
37  //! Assignment operator
38  GeneralMatrix& operator=(const GeneralMatrix& right);
39 
40  //! Destructor. Does nothing.
41  virtual ~GeneralMatrix();
42 
43  //! Duplicator member function
44  /*!
45  * This function will duplicate the matrix given a generic GeneralMatrix
46  *
47  * @return Returns a pointer to the malloced object
48  */
49  virtual GeneralMatrix* duplMyselfAsGeneralMatrix() const = 0;
50 
51  //! Zero the matrix elements
52  virtual void zero() = 0;
53 
54  //! Multiply A*b and write result to prod.
55  /*!
56  * @param b Vector to do the rh multiplication
57  * @param prod OUTPUT vector to receive the result
58  */
59  virtual void mult(const doublereal* b, doublereal* prod) const = 0;
60 
61  //! Multiply b*A and write result to prod.
62  /*!
63  * @param b Vector to do the lh multiplication
64  * @param prod OUTPUT vector to receive the result
65  */
66  virtual void leftMult(const doublereal* const b, doublereal* const prod) const = 0;
67 
68  //! Factors the A matrix, overwriting A.
69  /*
70  * We flip m_factored boolean to indicate that the matrix is now A-1.
71  */
72  virtual int factor() = 0;
73 
74  //! Factors the A matrix using the QR algorithm, overwriting A
75  /*!
76  * we set m_factored to 2 to indicate the matrix is now QR factored
77  *
78  * @return Returns the info variable from lapack
79  */
80  virtual int factorQR() = 0;
81 
82  //! Returns an estimate of the inverse of the condition number for the matrix
83  /*!
84  * The matrix must have been previously factored using the QR algorithm
85  *
86  * @return returns the inverse of the condition number
87  */
88  virtual doublereal rcondQR() = 0;
89 
90  //! Returns an estimate of the inverse of the condition number for the matrix
91  /*!
92  * The matrix must have been previously factored using the LU algorithm
93  *
94  * @param a1norm Norm of the matrix
95  *
96  * @return returns the inverse of the condition number
97  */
98  virtual doublereal rcond(doublereal a1norm) = 0;
99 
100  //! Change the way the matrix is factored
101  /*!
102  * @param fAlgorithm integer
103  * 0 LU factorization
104  * 1 QR factorization
105  */
106  virtual void useFactorAlgorithm(int fAlgorithm) = 0;
107 
108  //! Return the factor algorithm used
109  virtual int factorAlgorithm() const = 0;
110 
111  //! Calculate the one norm of the matrix
112  virtual doublereal oneNorm() const = 0;
113 
114  //! Return the number of rows in the matrix
115  virtual size_t nRows() const = 0;
116 
117  //! Return the size and structure of the matrix
118  /*!
119  * @param iStruct OUTPUT Pointer to a vector of ints that describe the structure of the matrix.
120  *
121  * @return returns the number of rows and columns in the matrix.
122  */
123  virtual size_t nRowsAndStruct(size_t* const iStruct = 0) const = 0;
124 
125  //! clear the factored flag
126  virtual void clearFactorFlag() = 0;
127 
128  //! Solves the Ax = b system returning x in the b spot.
129  /*!
130  * @param b Vector for the rhs of the equation system
131  */
132  virtual int solve(doublereal* b) = 0;
133 
134  //! true if the current factorization is up to date with the matrix
135  virtual bool factored() const = 0;
136 
137  //! Return a pointer to the top of column j, columns are assumed to be contiguous in memory
138  /*!
139  * @param j Value of the column
140  *
141  * @return Returns a pointer to the top of the column
142  */
143  virtual doublereal* ptrColumn(size_t j) = 0;
144 
145  //! Index into the (i,j) element
146  /*!
147  * @param i row
148  * @param j column
149  *
150  * Returns a changeable reference to the matrix entry
151  */
152  virtual doublereal& operator()(size_t i, size_t j) = 0;
153 
154  //! Constant Index into the (i,j) element
155  /*!
156  * @param i row
157  * @param j column
158  *
159  * Returns an unchangeable reference to the matrix entry
160  */
161  virtual doublereal operator()(size_t i, size_t j) const = 0;
162 
163  //! Copy the data from one array into another without doing any checking
164  /*!
165  * This differs from the assignment operator as no resizing is done and memcpy() is used.
166  * @param y Array to be copied
167  */
168  virtual void copyData(const GeneralMatrix& y) = 0;
169 
170  //! Return an iterator pointing to the first element
171  /*!
172  * We might drop this later
173  */
174  virtual vector_fp::iterator begin() = 0;
175 
176  //! Return a const iterator pointing to the first element
177  /*!
178  * We might drop this later
179  */
180  virtual vector_fp::const_iterator begin() const = 0;
181 
182  //! Return a vector of const pointers to the columns
183  /*!
184  * Note the value of the pointers are protected by their being const.
185  * However, the value of the matrix is open to being changed.
186  *
187  * @return returns a vector of pointers to the top of the columns
188  * of the matrices.
189  */
190  virtual doublereal* const* colPts() = 0;
191 
192  //! Check to see if we have any zero rows in the jacobian
193  /*!
194  * This utility routine checks to see if any rows are zero.
195  * The smallest row is returned along with the largest coefficient in that row
196  *
197  * @param valueSmall OUTPUT value of the largest coefficient in the smallest row
198  *
199  * @return index of the row that is most nearly zero
200  */
201  virtual size_t checkRows(doublereal& valueSmall) const = 0;
202 
203  //! Check to see if we have any zero columns in the jacobian
204  /*!
205  * This utility routine checks to see if any columns are zero.
206  * The smallest column is returned along with the largest coefficient in that column
207  *
208  * @param valueSmall OUTPUT value of the largest coefficient in the smallest column
209  *
210  * @return index of the column that is most nearly zero
211  */
212  virtual size_t checkColumns(doublereal& valueSmall) const = 0;
213 
214  //! Matrix type
215  /*!
216  * 0 Square
217  * 1 Banded
218  */
220 
221 };
222 }
223 #endif
virtual int factorQR()=0
Factors the A matrix using the QR algorithm, overwriting A.
virtual doublereal * ptrColumn(size_t j)=0
Return a pointer to the top of column j, columns are assumed to be contiguous in memory.
int matrixType_
Matrix type.
virtual void useFactorAlgorithm(int fAlgorithm)=0
Change the way the matrix is factored.
virtual GeneralMatrix * duplMyselfAsGeneralMatrix() const =0
Duplicator member function.
virtual doublereal oneNorm() const =0
Calculate the one norm of the matrix.
virtual size_t nRows() const =0
Return the number of rows in the matrix.
virtual doublereal rcond(doublereal a1norm)=0
Returns an estimate of the inverse of the condition number for the matrix.
GeneralMatrix(int matType)
Base Constructor.
virtual vector_fp::iterator begin()=0
Return an iterator pointing to the first element.
This file contains definitions of terms that are used in internal routines and are unlikely to need m...
Generic matrix.
Definition: GeneralMatrix.h:23
virtual int factorAlgorithm() const =0
Return the factor algorithm used.
virtual int factor()=0
Factors the A matrix, overwriting A.
virtual void zero()=0
Zero the matrix elements.
virtual doublereal & operator()(size_t i, size_t j)=0
Index into the (i,j) element.
virtual bool factored() const =0
true if the current factorization is up to date with the matrix
virtual ~GeneralMatrix()
Destructor. Does nothing.
virtual size_t checkRows(doublereal &valueSmall) const =0
Check to see if we have any zero rows in the jacobian.
virtual void clearFactorFlag()=0
clear the factored flag
virtual void leftMult(const doublereal *const b, doublereal *const prod) const =0
Multiply b*A and write result to prod.
virtual void mult(const doublereal *b, doublereal *prod) const =0
Multiply A*b and write result to prod.
virtual size_t nRowsAndStruct(size_t *const iStruct=0) const =0
Return the size and structure of the matrix.
virtual void copyData(const GeneralMatrix &y)=0
Copy the data from one array into another without doing any checking.
virtual size_t checkColumns(doublereal &valueSmall) const =0
Check to see if we have any zero columns in the jacobian.
virtual int solve(doublereal *b)=0
Solves the Ax = b system returning x in the b spot.
virtual doublereal *const * colPts()=0
Return a vector of const pointers to the columns.
GeneralMatrix & operator=(const GeneralMatrix &right)
Assignment operator.
virtual doublereal rcondQR()=0
Returns an estimate of the inverse of the condition number for the matrix.