Cantera  2.3.0
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 // This file is part of Cantera. See License.txt in the top-level directory or
8 // at http://www.cantera.org/license.txt for license and copyright information.
9 
10 #ifndef CT_GENERALMATRIX_H
11 #define CT_GENERALMATRIX_H
12 
13 #include "cantera/base/ct_defs.h"
15 #include "cantera/base/global.h"
16 
17 namespace Cantera
18 {
19 
20 //! Generic matrix
22 {
23 public:
24  //! Base Constructor
25  explicit GeneralMatrix(int matType=-1) : m_factored(false) {
26  if (matType != -1) {
27  warn_deprecated("GeneralMatrix", "Integer argument to constructor "
28  "is deprecated and will be removed after Cantera 2.3.");
29  }
30  }
31 
32  virtual ~GeneralMatrix() {}
33 
34  //! Duplicator member function
35  /*!
36  * This function will duplicate the matrix given a generic GeneralMatrix
37  *
38  * @returns a pointer to the newly created object
39  */
40  virtual GeneralMatrix* duplMyselfAsGeneralMatrix() const = 0;
41 
42  //! Zero the matrix elements
43  virtual void zero() = 0;
44 
45  //! Multiply A*b and write result to prod.
46  /*!
47  * @param b Vector to do the rh multiplication
48  * @param prod OUTPUT vector to receive the result
49  */
50  virtual void mult(const doublereal* b, doublereal* prod) const = 0;
51 
52  //! Multiply b*A and write result to prod.
53  /*!
54  * @param b Vector to do the lh multiplication
55  * @param prod OUTPUT vector to receive the result
56  */
57  virtual void leftMult(const doublereal* const b, doublereal* const prod) const = 0;
58 
59  //! Factors the A matrix, overwriting A.
60  /*!
61  * We flip m_factored boolean to indicate that the matrix is now A-1.
62  */
63  virtual int factor() = 0;
64 
65  //! Factors the A matrix using the QR algorithm, overwriting A
66  /*!
67  * we set m_factored to 2 to indicate the matrix is now QR factored
68  *
69  * @returns the info variable from LAPACK
70  */
71  virtual int factorQR() {
72  throw NotImplementedError("GeneralMatrix::factorQR");
73  }
74 
75  //! Returns an estimate of the inverse of the condition number for the matrix
76  /*!
77  * The matrix must have been previously factored using the QR algorithm
78  *
79  * @returns the inverse of the condition number
80  */
81  virtual doublereal rcondQR() {
82  throw NotImplementedError("GeneralMatrix::rcondQR");
83  }
84 
85  //! Returns an estimate of the inverse of the condition number for the matrix
86  /*!
87  * The matrix must have been previously factored using the LU algorithm
88  *
89  * @param a1norm Norm of the matrix
90  * @returns the inverse of the condition number
91  */
92  virtual doublereal rcond(doublereal a1norm) = 0;
93 
94  //! Change the way the matrix is factored
95  /*!
96  * @param fAlgorithm integer
97  * 0 LU factorization
98  * 1 QR factorization
99  */
100  virtual void useFactorAlgorithm(int fAlgorithm) {
101  throw NotImplementedError("GeneralMatrix::useFactorAlgorithm");
102  };
103 
104  //! Return the factor algorithm used
105  virtual int factorAlgorithm() const = 0;
106 
107  //! Calculate the one norm of the matrix
108  virtual doublereal oneNorm() const = 0;
109 
110  //! Return the number of rows in the matrix
111  virtual size_t nRows() const = 0;
112 
113  //! Return the size and structure of the matrix
114  /*!
115  * @param iStruct OUTPUT Pointer to a vector of ints that describe the
116  * structure of the matrix.
117  * @returns the number of rows and columns in the matrix.
118  * @deprecated Unused. To be removed after Cantera 2.3.
119  */
120  virtual size_t nRowsAndStruct(size_t* const iStruct = 0) const = 0;
121 
122  //! clear the factored flag
123  virtual void clearFactorFlag() {
124  m_factored = 0;
125  };
126 
127  //! Solves the Ax = b system returning x in the b spot.
128  /*!
129  * @param b Vector for the RHS of the equation system
130  * @param nrhs Number of right-hand sides to solve, default 1
131  * @param ldb Leading dimension of the right-hand side array. Defaults to
132  * nRows()
133  */
134  virtual int solve(doublereal* b, size_t nrhs=1, size_t ldb=0) = 0;
135 
136  //! true if the current factorization is up to date with the matrix
137  virtual bool factored() const {
138  return (m_factored != 0);
139  }
140 
141  //! Return a pointer to the top of column j, columns are assumed to be
142  //! contiguous in memory
143  /*!
144  * @param j Value of the column
145  * @returns a pointer to the top of the column
146  */
147  virtual doublereal* ptrColumn(size_t j) = 0;
148 
149  //! Index into the (i,j) element
150  /*!
151  * @param i row
152  * @param j column
153  * @returns a changeable reference to the matrix entry
154  */
155  virtual doublereal& operator()(size_t i, size_t j) = 0;
156 
157  //! Constant Index into the (i,j) element
158  /*!
159  * @param i row
160  * @param j column
161  * @returns an unchangeable reference to the matrix entry
162  */
163  virtual doublereal operator()(size_t i, size_t j) const = 0;
164 
165  //! Return an iterator pointing to the first element
166  /*!
167  * We might drop this later
168  */
169  virtual vector_fp::iterator begin() = 0;
170 
171  //! Return a const iterator pointing to the first element
172  /*!
173  * We might drop this later
174  */
175  virtual vector_fp::const_iterator begin() const = 0;
176 
177  //! Return a vector of const pointers to the columns
178  /*!
179  * Note the value of the pointers are protected by their being const.
180  * However, the value of the matrix is open to being changed.
181  *
182  * @returns a vector of pointers to the top of the columns of the matrices.
183  */
184  virtual doublereal* const* colPts() = 0;
185 
186  //! Check to see if we have any zero rows in the Jacobian
187  /*!
188  * This utility routine checks to see if any rows are zero. The smallest row
189  * is returned along with the largest coefficient in that row
190  *
191  * @param valueSmall OUTPUT value of the largest coefficient in the smallest row
192  * @return index of the row that is most nearly zero
193  */
194  virtual size_t checkRows(doublereal& valueSmall) const = 0;
195 
196  //! Check to see if we have any zero columns in the Jacobian
197  /*!
198  * This utility routine checks to see if any columns are zero. The smallest
199  * column is returned along with the largest coefficient in that column
200  *
201  * @param valueSmall OUTPUT value of the largest coefficient in the smallest column
202  * @return index of the column that is most nearly zero
203  */
204  virtual size_t checkColumns(doublereal& valueSmall) const = 0;
205 
206 protected:
207  //! Indicates whether the matrix is factored. 0 for unfactored; Non-zero
208  //! values indicate a particular factorization (LU=1, QR=2).
210 };
211 
212 }
213 #endif
virtual doublereal * ptrColumn(size_t j)=0
Return a pointer to the top of column j, columns are assumed to be contiguous in memory.
An error indicating that an unimplemented function has been called.
Definition: ctexceptions.h:193
virtual GeneralMatrix * duplMyselfAsGeneralMatrix() const =0
Duplicator member function.
int m_factored
Indicates whether the matrix is factored.
virtual doublereal oneNorm() const =0
Calculate the one norm of the matrix.
virtual void clearFactorFlag()
clear the factored flag
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.
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...
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
Definition: global.cpp:54
Generic matrix.
Definition: GeneralMatrix.h:21
virtual int factorQR()
Factors the A matrix using the QR algorithm, overwriting A.
Definition: GeneralMatrix.h:71
virtual int factorAlgorithm() const =0
Return the factor algorithm used.
This file contains definitions for utility functions and text for modules, inputfiles, logs, textlogs, (see Input File Handling, Diagnostic Output, and Writing messages to the screen).
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
true if the current factorization is up to date with the matrix
virtual size_t checkRows(doublereal &valueSmall) const =0
Check to see if we have any zero rows in the Jacobian.
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.
GeneralMatrix(int matType=-1)
Base Constructor.
Definition: GeneralMatrix.h:25
virtual int solve(doublereal *b, size_t nrhs=1, size_t ldb=0)=0
Solves the Ax = b system returning x in the b spot.
virtual size_t nRowsAndStruct(size_t *const iStruct=0) const =0
Return the size and structure of the matrix.
virtual size_t checkColumns(doublereal &valueSmall) const =0
Check to see if we have any zero columns in the Jacobian.
virtual doublereal *const * colPts()=0
Return a vector of const pointers to the columns.
Namespace for the Cantera kernel.
Definition: application.cpp:29
virtual doublereal rcondQR()
Returns an estimate of the inverse of the condition number for the matrix.
Definition: GeneralMatrix.h:81
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
virtual void useFactorAlgorithm(int fAlgorithm)
Change the way the matrix is factored.