Cantera  2.1.2
BandMatrix.h
Go to the documentation of this file.
1 /**
2  * @file BandMatrix.h
3  * Declarations for the class BandMatrix
4  * which is a child class of GeneralMatrix for banded matrices handled by solvers
5  * (see class \ref numerics and \link Cantera::BandMatrix BandMatrix\endlink).
6  */
7 
8 // Copyright 2001 California Institute of Technology
9 
10 #ifndef CT_BANDMATRIX_H
11 #define CT_BANDMATRIX_H
12 
13 #include "cantera/base/ct_defs.h"
14 #include "ctlapack.h"
15 #include "cantera/base/utilities.h"
17 #include "GeneralMatrix.h"
18 
19 namespace Cantera
20 {
21 
22 //! A class for banded matrices, involving matrix inversion processes.
23 //! The class is based upon the LAPACK banded storage matrix format.
24 /*!
25  * An important issue with this class is that it stores both the original data
26  * and the LU factorization of the data. This means that the banded matrix typically
27  * will take up twice the room that it is expected to take.
28  *
29  * QR factorizations of banded matrices are not included in the original LAPACK work.
30  * Add-ons are available. However, they are not included here. Instead we just use the
31  * stock LU decompositions.
32  *
33  * This class is a derived class of the base class GeneralMatrix. However, within
34  * the oneD directory, the class is used as is, without reference to the GeneralMatrix
35  * base type.
36  */
37 class BandMatrix : public GeneralMatrix
38 {
39 
40 public:
41 
42  //! Base Constructor
43  /*!
44  * * Create an \c 0 by \c 0 matrix, and initialize all elements to \c 0.
45  */
46  BandMatrix();
47 
48  //! Creates a banded matrix and sets all elements to zero
49  /*!
50  * Create an \c n by \c n banded matrix, and initialize all elements to \c v.
51  *
52  * @param n size of the square matrix
53  * @param kl band size on the lower portion of the matrix
54  * @param ku band size on the upper portion of the matrix
55  * @param v initial value of all matrix components.
56  */
57  BandMatrix(size_t n, size_t kl, size_t ku, doublereal v = 0.0);
58 
59  //! Copy constructor
60  /*!
61  * @param y Matrix to be copied
62  */
63  BandMatrix(const BandMatrix& y);
64 
65  //! assignment operator
66  /*!
67  * @param y reference to the matrix to be copied
68  */
69  BandMatrix& operator=(const BandMatrix& y);
70 
71  //! Resize the matrix problem
72  /*!
73  * All data is lost
74  *
75  * @param n size of the square matrix
76  * @param kl band size on the lower portion of the matrix
77  * @param ku band size on the upper portion of the matrix
78  * @param v initial value of all matrix components.
79  */
80  void resize(size_t n, size_t kl, size_t ku, doublereal v = 0.0);
81 
82  //! Fill or zero the matrix
83  /*!
84  * @param v Fill value, defaults to zero.
85  */
86  void bfill(doublereal v = 0.0);
87 
88  doublereal& operator()(size_t i, size_t j);
89  doublereal operator()(size_t i, size_t j) const;
90 
91  //! Return a changeable reference to element (i,j).
92  /*!
93  * Since this method may alter the element value, it may need to be refactored, so
94  * the flag m_factored is set to false.
95  *
96  * @param i row
97  * @param j column
98  *
99  * @return Returns a reference to the value of the matrix entry
100  */
101  doublereal& value(size_t i, size_t j);
102 
103  //! Return the value of element (i,j).
104  /*!
105  * This method does not alter the array.
106  * @param i row
107  * @param j column
108  *
109  * @return Returns the value of the matrix entry
110  */
111  doublereal value(size_t i, size_t j) const;
112 
113  //! Returns the location in the internal 1D array corresponding to the (i,j) element in the banded array
114  /*!
115  * @param i row
116  * @param j column
117  *
118  * @return Returns the index of the matrix entry
119  */
120  size_t index(size_t i, size_t j) const;
121 
122  //! Return the value of the (i,j) element for (i,j) within the bandwidth.
123  /*!
124  * For efficiency, this method does not check that (i,j) are within the bandwidth; it is up to the calling
125  * program to insure that this is true.
126  *
127  * @param i row
128  * @param j column
129  *
130  * @return Returns the value of the matrix entry
131  */
132  doublereal _value(size_t i, size_t j) const;
133 
134  virtual size_t nRows() const;
135 
136  //! Return the size and structure of the matrix
137  /*!
138  * @param iStruct OUTPUT Pointer to a vector of ints that describe the structure of the matrix.
139  * istruct[0] = kl
140  * istruct[1] = ku
141  *
142  * @return returns the number of rows and columns in the matrix.
143  */
144  virtual size_t nRowsAndStruct(size_t* const iStruct = 0) const;
145 
146  //! Number of columns
147  size_t nColumns() const;
148 
149  //! Number of subdiagonals
150  size_t nSubDiagonals() const;
151 
152  //! Number of superdiagonals
153  size_t nSuperDiagonals() const;
154 
155  //! Return the number of rows of storage needed for the band storage
156  size_t ldim() const;
157 
158  //! Return a reference to the pivot vector
159  vector_int& ipiv();
160 
161  virtual void mult(const doublereal* b, doublereal* prod) const;
162  virtual void leftMult(const doublereal* const b, doublereal* const prod) const;
163 
164  //! Perform an LU decomposition, the LAPACK routine DGBTRF is used.
165  /*!
166  * The factorization is saved in ludata.
167  *
168  * @return Return a success flag.
169  * 0 indicates a success
170  * ~0 Some error occurred, see the LAPACK documentation
171  */
172  int factor();
173 
174  //! Solve the matrix problem Ax = b
175  /*!
176  * @param b INPUT rhs of the problem
177  * @param x OUTPUT solution to the problem
178  *
179  * @return Return a success flag
180  * 0 indicates a success
181  * ~0 Some error occurred, see the LAPACK documentation
182  */
183  int solve(const doublereal* const b, doublereal* const x);
184 
185  //! Solve the matrix problem Ax = b
186  /*!
187  * @param b INPUT rhs of the problem
188  * OUTPUT solution to the problem
189  *
190  * @return Return a success flag
191  * 0 indicates a success
192  * ~0 Some error occurred, see the LAPACK documentation
193  */
194  int solve(doublereal* b);
195 
196  //! Returns an iterator for the start of the band storage data
197  /*!
198  * Iterator points to the beginning of the data, and it is changeable.
199  */
200  virtual vector_fp::iterator begin();
201 
202  //! Returns an iterator for the end of the band storage data
203  /*!
204  * Iterator points to the end of the data, and it is changeable.
205  */
206  vector_fp::iterator end();
207 
208  //! Returns a const iterator for the start of the band storage data
209  /*!
210  * Iterator points to the beginning of the data, and it is not changeable.
211  */
212  vector_fp::const_iterator begin() const;
213 
214  //! Returns a const iterator for the end of the band storage data
215  /*!
216  * Iterator points to the end of the data, and it is not changeable.
217  */
218  vector_fp::const_iterator end() const;
219 
220  virtual void zero();
221 
222  //! Factors the A matrix using the QR algorithm, overwriting A
223  /*!
224  * we set m_factored to 2 to indicate the matrix is now QR factored
225  *
226  * @return Returns the info variable from lapack
227  */
228  virtual int factorQR();
229 
230  //! Returns an estimate of the inverse of the condition number for the matrix
231  /*!
232  * The matrix must have been previously factored using the QR algorithm
233  *
234  * @return returns the inverse of the condition number
235  */
236  virtual doublereal rcondQR();
237 
238  //! Returns an estimate of the inverse of the condition number for the matrix
239  /*!
240  * The matrix must have been previously factored using the LU algorithm
241  *
242  * @param a1norm Norm of the matrix
243  *
244  * @return returns the inverse of the condition number
245  */
246  virtual doublereal rcond(doublereal a1norm);
247 
248  //! Change the way the matrix is factored
249  /*!
250  * @param fAlgorithm integer
251  * 0 LU factorization
252  * 1 QR factorization
253  */
254  virtual void useFactorAlgorithm(int fAlgorithm);
255 
256  //! Returns the factor algorithm used
257  /*!
258  * 0 LU decomposition
259  * 1 QR decomposition
260  *
261  * This routine will always return 0
262  */
263  virtual int factorAlgorithm() const;
264 
265  //! Returns the one norm of the matrix
266  virtual doublereal oneNorm() const;
267 
268  virtual GeneralMatrix* duplMyselfAsGeneralMatrix() const;
269 
270  //! Report whether the current matrix has been factored.
271  virtual bool factored() const;
272 
273  //! Return a pointer to the top of column j, column values are assumed to be contiguous in memory
274  /*!
275  * The LAPACK bandstructure has column values which are contiguous in memory:
276  *
277  * On entry, the matrix A in band storage, in rows KL+1 to
278  * 2*KL+KU+1; rows 1 to KL of the array need not be set.
279  * The j-th column of A is stored in the j-th column of the
280  * array AB as follows:
281  * AB(KL + KU + 1 + i - j,j) = A(i,j) for max(1, j - KU) <= i <= min(m, j + KL)
282  *
283  * This routine returns the position of AB(1,j) (fortran-1 indexing) in the above format
284  *
285  * So to address the (i,j) position, you use the following indexing:
286  *
287  * double *colP_j = matrix.ptrColumn(j);
288  * double a_i_j = colP_j[kl + ku + i - j];
289  *
290  * @param j Value of the column
291  *
292  * @return Returns a pointer to the top of the column
293  */
294  virtual doublereal* ptrColumn(size_t j);
295 
296  //! Return a vector of const pointers to the columns
297  /*!
298  * Note the value of the pointers are protected by their being const.
299  * However, the value of the matrix is open to being changed.
300  *
301  * @return returns a vector of pointers to the top of the columns
302  * of the matrices.
303  */
304  virtual doublereal* const* colPts();
305 
306  //! Copy the data from one array into another without doing any checking
307  /*!
308  * This differs from the assignment operator as no resizing is done and memcpy() is used.
309  * @param y Array to be copied
310  */
311  virtual void copyData(const GeneralMatrix& y);
312 
313  //! Clear the factored flag
314  virtual void clearFactorFlag();
315 
316  //! Check to see if we have any zero rows in the jacobian
317  /*!
318  * This utility routine checks to see if any rows are zero.
319  * The smallest row is returned along with the largest coefficient in that row
320  *
321  * @param valueSmall OUTPUT value of the largest coefficient in the smallest row
322  *
323  * @return index of the row that is most nearly zero
324  */
325  virtual size_t checkRows(doublereal& valueSmall) const;
326 
327  //! Check to see if we have any zero columns in the jacobian
328  /*!
329  * This utility routine checks to see if any columns are zero.
330  * The smallest column is returned along with the largest coefficient in that column
331  *
332  * @param valueSmall OUTPUT value of the largest coefficient in the smallest column
333  *
334  * @return index of the column that is most nearly zero
335  */
336  virtual size_t checkColumns(doublereal& valueSmall) const;
337 
338 protected:
339 
340  //! Matrix data
342 
343  //! Factorized data
345 
346  //! Boolean indicating whether the matrix is factored
348 
349  //! Number of rows and columns of the matrix
350  size_t m_n;
351 
352  //! Number of subdiagonals of the matrix
353  size_t m_kl;
354 
355  //! Number of super diagonals of the matrix
356  size_t m_ku;
357 
358  //! value of zero
359  doublereal m_zero;
360 
361  //! Pivot vector
363 
364  //! Vector of column pointers
365  std::vector<doublereal*> m_colPtrs;
366 
367  //! Extra work array needed - size = n
369 
370  //! Extra dp work array needed - size = 3n
372 
373 private:
374 
375  //! Error function that gets called for unhandled cases
376  /*!
377  * @param msg String containing the message.
378  */
379  void err(const std::string& msg) const;
380 };
381 
382 //! Utility routine to print out the matrix
383 /*!
384  * @param s ostream to print the matrix out to
385  * @param m Matrix to be printed
386  *
387  * @return Returns a reference to the ostream
388  */
389 std::ostream& operator<<(std::ostream& s, const BandMatrix& m);
390 
391 }
392 
393 #endif
size_t m_kl
Number of subdiagonals of the matrix.
Definition: BandMatrix.h:353
virtual vector_fp::iterator begin()
Returns an iterator for the start of the band storage data.
Definition: BandMatrix.cpp:293
vector_int iwork_
Extra work array needed - size = n.
Definition: BandMatrix.h:368
virtual doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j, column values are assumed to be contiguous in memory...
Definition: BandMatrix.cpp:466
vector_fp work_
Extra dp work array needed - size = 3n.
Definition: BandMatrix.h:371
virtual int factorQR()
Factors the A matrix using the QR algorithm, overwriting A.
Definition: BandMatrix.cpp:333
BandMatrix()
Base Constructor.
Definition: BandMatrix.cpp:24
vector_int & ipiv()
Return a reference to the pivot vector.
Definition: BandMatrix.cpp:205
BandMatrix & operator=(const BandMatrix &y)
assignment operator
Definition: BandMatrix.cpp:78
doublereal _value(size_t i, size_t j) const
Return the value of the (i,j) element for (i,j) within the bandwidth.
Definition: BandMatrix.cpp:163
Various templated functions that carry out common vector operations (see Templated Utility Functions)...
size_t nSuperDiagonals() const
Number of superdiagonals.
Definition: BandMatrix.cpp:195
virtual size_t checkRows(doublereal &valueSmall) const
Check to see if we have any zero rows in the jacobian.
Definition: BandMatrix.cpp:404
int solve(const doublereal *const b, doublereal *const x)
Solve the matrix problem Ax = b.
Definition: BandMatrix.cpp:267
doublereal & operator()(size_t i, size_t j)
Index into the (i,j) element.
Definition: BandMatrix.cpp:128
This file contains definitions of terms that are used in internal routines and are unlikely to need m...
Generic matrix.
Definition: GeneralMatrix.h:23
size_t nColumns() const
Number of columns.
Definition: BandMatrix.cpp:183
vector_fp data
Matrix data.
Definition: BandMatrix.h:341
size_t m_n
Number of rows and columns of the matrix.
Definition: BandMatrix.h:350
void resize(size_t n, size_t kl, size_t ku, doublereal v=0.0)
Resize the matrix problem.
Definition: BandMatrix.cpp:99
virtual void leftMult(const doublereal *const b, doublereal *const prod) const
Multiply b*A and write result to prod.
Definition: BandMatrix.cpp:230
virtual GeneralMatrix * duplMyselfAsGeneralMatrix() const
Duplicator member function.
Definition: BandMatrix.cpp:456
virtual size_t checkColumns(doublereal &valueSmall) const
Check to see if we have any zero columns in the jacobian.
Definition: BandMatrix.cpp:430
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:167
virtual void clearFactorFlag()
Clear the factored flag.
Definition: BandMatrix.cpp:484
virtual size_t nRowsAndStruct(size_t *const iStruct=0) const
Return the size and structure of the matrix.
Definition: BandMatrix.cpp:173
virtual doublereal rcondQR()
Returns an estimate of the inverse of the condition number for the matrix.
Definition: BandMatrix.cpp:339
void err(const std::string &msg) const
Error function that gets called for unhandled cases.
Definition: BandMatrix.cpp:328
virtual doublereal rcond(doublereal a1norm)
Returns an estimate of the inverse of the condition number for the matrix.
Definition: BandMatrix.cpp:345
vector_int m_ipiv
Pivot vector.
Definition: BandMatrix.h:362
bool m_factored
Boolean indicating whether the matrix is factored.
Definition: BandMatrix.h:347
doublereal & value(size_t i, size_t j)
Return a changeable reference to element (i,j).
Definition: BandMatrix.cpp:138
size_t ldim() const
Return the number of rows of storage needed for the band storage.
Definition: BandMatrix.cpp:200
virtual doublereal *const * colPts()
Return a vector of const pointers to the columns.
Definition: BandMatrix.cpp:471
virtual void copyData(const GeneralMatrix &y)
Copy the data from one array into another without doing any checking.
Definition: BandMatrix.cpp:476
virtual size_t nRows() const
Return the number of rows in the matrix.
Definition: BandMatrix.cpp:168
vector_fp ludata
Factorized data.
Definition: BandMatrix.h:344
virtual int factorAlgorithm() const
Returns the factor algorithm used.
Definition: BandMatrix.cpp:381
size_t index(size_t i, size_t j) const
Returns the location in the internal 1D array corresponding to the (i,j) element in the banded array...
Definition: BandMatrix.cpp:155
size_t m_ku
Number of super diagonals of the matrix.
Definition: BandMatrix.h:356
virtual bool factored() const
Report whether the current matrix has been factored.
Definition: BandMatrix.cpp:461
virtual void useFactorAlgorithm(int fAlgorithm)
Change the way the matrix is factored.
Definition: BandMatrix.cpp:376
void bfill(doublereal v=0.0)
Fill or zero the matrix.
Definition: BandMatrix.cpp:116
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
doublereal m_zero
value of zero
Definition: BandMatrix.h:359
virtual doublereal oneNorm() const
Returns the one norm of the matrix.
Definition: BandMatrix.cpp:386
int factor()
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
Definition: BandMatrix.cpp:248
virtual void mult(const doublereal *b, doublereal *prod) const
Multiply A*b and write result to prod.
Definition: BandMatrix.cpp:213
Declarations for the class GeneralMatrix which is a virtual base class for matrices handled by solver...
virtual void zero()
Zero the matrix elements.
Definition: BandMatrix.cpp:122
std::vector< doublereal * > m_colPtrs
Vector of column pointers.
Definition: BandMatrix.h:365
size_t nSubDiagonals() const
Number of subdiagonals.
Definition: BandMatrix.cpp:189
Definitions for the classes that are thrown when Cantera experiences an error condition (also contain...
vector_fp::iterator end()
Returns an iterator for the end of the band storage data.
Definition: BandMatrix.cpp:299
A class for banded matrices, involving matrix inversion processes.
Definition: BandMatrix.h:37