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