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