Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BandMatrix.cpp
Go to the documentation of this file.
1 /**
2  * @file BandMatrix.cpp
3  *
4  * Banded matrices.
5  */
6 
7 // Copyright 2001 California Institute of Technology
8 
11 #include "cantera/base/utilities.h"
13 
14 #include <cstring>
15 #include <fstream>
16 
17 using namespace std;
18 
19 namespace Cantera
20 {
21 
22 BandMatrix::BandMatrix() :
23  GeneralMatrix(1),
24  m_n(0),
25  m_kl(0),
26  m_ku(0),
27  m_zero(0.0)
28 {
29 }
30 
31 BandMatrix::BandMatrix(size_t n, size_t kl, size_t ku, doublereal v) :
32  GeneralMatrix(1),
33  m_n(n),
34  m_kl(kl),
35  m_ku(ku),
36  m_zero(0.0)
37 {
38  data.resize(n*(2*kl + ku + 1));
39  ludata.resize(n*(2*kl + ku + 1));
40  fill(data.begin(), data.end(), v);
41  fill(ludata.begin(), ludata.end(), 0.0);
42  m_ipiv.resize(m_n);
43  m_colPtrs.resize(n);
44  size_t ldab = (2*kl + ku + 1);
45  for (size_t j = 0; j < n; j++) {
46  m_colPtrs[j] = &(data[ldab * j]);
47  }
48 }
49 
51  GeneralMatrix(y),
52  m_n(0),
53  m_kl(0),
54  m_ku(0),
55  m_zero(0.0)
56 {
57  m_n = y.m_n;
58  m_kl = y.m_kl;
59  m_ku = y.m_ku;
60  data = y.data;
61  ludata = y.ludata;
62  m_ipiv = y.m_ipiv;
63  m_colPtrs.resize(m_n);
64  size_t ldab = (2 *m_kl + m_ku + 1);
65  for (size_t j = 0; j < m_n; j++) {
66  m_colPtrs[j] = &(data[ldab * j]);
67  }
68 }
69 
71 {
72  if (&y == this) {
73  return *this;
74  }
76  m_n = y.m_n;
77  m_kl = y.m_kl;
78  m_ku = y.m_ku;
79  m_ipiv = y.m_ipiv;
80  data = y.data;
81  ludata = y.ludata;
82  m_colPtrs.resize(m_n);
83  size_t ldab = (2 * m_kl + m_ku + 1);
84  for (size_t j = 0; j < m_n; j++) {
85  m_colPtrs[j] = &(data[ldab * j]);
86  }
87  return *this;
88 }
89 
90 void BandMatrix::resize(size_t n, size_t kl, size_t ku, doublereal v)
91 {
92  m_n = n;
93  m_kl = kl;
94  m_ku = ku;
95  data.resize(n*(2*kl + ku + 1));
96  ludata.resize(n*(2*kl + ku + 1));
97  m_ipiv.resize(m_n);
98  fill(data.begin(), data.end(), v);
99  m_colPtrs.resize(m_n);
100  size_t ldab = (2 * m_kl + m_ku + 1);
101  for (size_t j = 0; j < n; j++) {
102  m_colPtrs[j] = &(data[ldab * j]);
103  }
104  m_factored = false;
105 }
106 
107 void BandMatrix::bfill(doublereal v)
108 {
109  std::fill(data.begin(), data.end(), v);
110  m_factored = false;
111 }
112 
114 {
115  std::fill(data.begin(), data.end(), 0.0);
116  m_factored = false;
117 }
118 
119 doublereal& BandMatrix::operator()(size_t i, size_t j)
120 {
121  return value(i,j);
122 }
123 
124 doublereal BandMatrix::operator()(size_t i, size_t j) const
125 {
126  return value(i,j);
127 }
128 
129 doublereal& BandMatrix::value(size_t i, size_t j)
130 {
131  m_factored = false;
132  if (i + m_ku < j || i > j + m_kl) {
133  return m_zero;
134  }
135  return data[index(i,j)];
136 }
137 
138 doublereal BandMatrix::value(size_t i, size_t j) const
139 {
140  if (i + m_ku < j || i > j + m_kl) {
141  return 0.0;
142  }
143  return data[index(i,j)];
144 }
145 
146 size_t BandMatrix::index(size_t i, size_t j) const
147 {
148  int jj = static_cast<int>(j);
149  int ii = static_cast<int>(i);
150  size_t rw = (int) m_kl + (int) m_ku + (int) ii - jj;
151  return (2*m_kl + m_ku + 1)*j + rw;
152 }
153 
154 doublereal BandMatrix::_value(size_t i, size_t j) const
155 {
156  return data[index(i,j)];
157 }
158 
159 size_t BandMatrix::nRows() const
160 {
161  return m_n;
162 }
163 
164 size_t BandMatrix::nRowsAndStruct(size_t* const iStruct) const
165 {
166  if (iStruct) {
167  iStruct[0] = m_kl;
168  iStruct[1] = m_ku;
169  }
170  return m_n;
171 }
172 
173 size_t BandMatrix::nColumns() const
174 {
175  return m_n;
176 }
177 
179 {
180  return m_kl;
181 }
182 
184 {
185  return m_ku;
186 }
187 
188 size_t BandMatrix::ldim() const
189 {
190  return 2*m_kl + m_ku + 1;
191 }
192 
194 {
195  return m_ipiv;
196 }
197 
198 void BandMatrix::mult(const doublereal* b, doublereal* prod) const
199 {
200  int kl = static_cast<int>(m_kl);
201  int ku = static_cast<int>(m_ku);
202  int nr = static_cast<int>(nRows());
203  doublereal sum = 0.0;
204  for (int m = 0; m < nr; m++) {
205  sum = 0.0;
206  for (int j = m - kl; j <= m + ku; j++) {
207  if (j >= 0 && j < (int) m_n) {
208  sum += _value(m,j) * b[j];
209  }
210  }
211  prod[m] = sum;
212  }
213 }
214 
215 void BandMatrix::leftMult(const doublereal* const b, doublereal* const prod) const
216 {
217  int kl = static_cast<int>(m_kl);
218  int ku = static_cast<int>(m_ku);
219  int nc = static_cast<int>(nColumns());
220  doublereal sum = 0.0;
221  for (int n = 0; n < nc; n++) {
222  sum = 0.0;
223  for (int i = n - ku; i <= n + kl; i++) {
224  if (i >= 0 && i < (int) m_n) {
225  size_t ii = i;
226  sum += _value(ii,n) * b[ii];
227  }
228  }
229  prod[n] = sum;
230  }
231 }
232 
234 {
235  int info=0;
236  copy(data.begin(), data.end(), ludata.begin());
237  ct_dgbtrf(nRows(), nColumns(), nSubDiagonals(), nSuperDiagonals(),
238  DATA_PTR(ludata), ldim(), DATA_PTR(ipiv()), info);
239 
240  // if info = 0, LU decomp succeeded.
241  if (info == 0) {
242  m_factored = true;
243  } else {
244  m_factored = false;
245  ofstream fout("bandmatrix.csv");
246  fout << *this << endl;
247  fout.close();
248  }
249  return info;
250 }
251 
252 int BandMatrix::solve(const doublereal* const b, doublereal* const x)
253 {
254  copy(b, b + m_n, x);
255  return solve(x);
256 }
257 
258 int BandMatrix::solve(doublereal* b, size_t nrhs, size_t ldb)
259 {
260  int info = 0;
261  if (!m_factored) {
262  info = factor();
263  }
264  if (ldb == 0) {
265  ldb = nColumns();
266  }
267  if (info == 0)
268  ct_dgbtrs(ctlapack::NoTranspose, nColumns(), nSubDiagonals(),
269  nSuperDiagonals(), nrhs, DATA_PTR(ludata), ldim(),
270  DATA_PTR(ipiv()), b, ldb, info);
271 
272  // error handling
273  if (info != 0) {
274  ofstream fout("bandmatrix.csv");
275  fout << *this << endl;
276  fout.close();
277  }
278  return info;
279 }
280 
281 vector_fp::iterator BandMatrix::begin()
282 {
283  m_factored = false;
284  return data.begin();
285 }
286 
287 vector_fp::iterator BandMatrix::end()
288 {
289  m_factored = false;
290  return data.end();
291 }
292 
293 vector_fp::const_iterator BandMatrix::begin() const
294 {
295  return data.begin();
296 }
297 
298 vector_fp::const_iterator BandMatrix::end() const
299 {
300  return data.end();
301 }
302 
303 ostream& operator<<(ostream& s, const BandMatrix& m)
304 {
305  size_t nr = m.nRows();
306  size_t nc = m.nColumns();
307  for (size_t i = 0; i < nr; i++) {
308  for (size_t j = 0; j < nc; j++) {
309  s << m(i,j) << ", ";
310  }
311  s << endl;
312  }
313  return s;
314 }
315 
316 doublereal BandMatrix::rcond(doublereal a1norm)
317 {
318  int printLevel = 0;
319  int useReturnErrorCode = 0;
320  if (iwork_.size() < m_n) {
321  iwork_.resize(m_n);
322  }
323  if (work_.size() < 3 * m_n) {
324  work_.resize(3 * m_n);
325  }
326  doublereal rcond = 0.0;
327  if (m_factored != 1) {
328  throw CanteraError("BandMatrix::rcond()", "matrix isn't factored correctly");
329  }
330 
331  size_t ldab = (2 *m_kl + m_ku + 1);
332  int rinfo = 0;
333  rcond = ct_dgbcon('1', m_n, m_kl, m_ku, DATA_PTR(ludata), ldab, DATA_PTR(m_ipiv), a1norm, DATA_PTR(work_),
334  DATA_PTR(iwork_), rinfo);
335  if (rinfo != 0) {
336  if (printLevel) {
337  writelogf("BandMatrix::rcond(): DGBCON returned INFO = %d\n", rinfo);
338  }
339  if (! useReturnErrorCode) {
340  throw CanteraError("BandMatrix::rcond()", "DGBCON returned INFO = " + int2str(rinfo));
341  }
342  }
343  return rcond;
344 }
345 
347 {
348  return 0;
349 }
350 
351 doublereal BandMatrix::oneNorm() const
352 {
353  int ku = static_cast<int>(m_ku);
354  int kl = static_cast<int>(m_kl);
355  doublereal value = 0.0;
356  for (int j = 0; j < (int) m_n; j++) {
357  doublereal sum = 0.0;
358  doublereal* colP = m_colPtrs[j];
359  for (int i = j - ku; i <= j + kl; i++) {
360  sum += fabs(colP[kl + ku + i - j]);
361  }
362  value = std::max(sum, value);
363  }
364  return value;
365 }
366 
367 size_t BandMatrix::checkRows(doublereal& valueSmall) const
368 {
369  valueSmall = 1.0E300;
370  size_t iSmall = npos;
371  for (int i = 0; i < (int) m_n; i++) {
372  double valueS = 0.0;
373  for (int j = i - (int) m_kl; j <= i + (int) m_ku; j++) {
374  if (j >= 0 && j < (int) m_n) {
375  valueS = std::max(fabs(value(i,j)), valueS);
376  }
377  }
378  if (valueS < valueSmall) {
379  iSmall = i;
380  valueSmall = valueS;
381  if (valueSmall == 0.0) {
382  return iSmall;
383  }
384  }
385  }
386  return iSmall;
387 }
388 
389 size_t BandMatrix::checkColumns(doublereal& valueSmall) const
390 {
391  valueSmall = 1.0E300;
392  size_t jSmall = npos;
393  for (int j = 0; j < (int) m_n; j++) {
394  double valueS = 0.0;
395  for (int i = j - (int) m_ku; i <= j + (int) m_kl; i++) {
396  if (i >= 0 && i < (int) m_n) {
397  valueS = std::max(fabs(value(i,j)), valueS);
398  }
399  }
400  if (valueS < valueSmall) {
401  jSmall = j;
402  valueSmall = valueS;
403  if (valueSmall == 0.0) {
404  return jSmall;
405  }
406  }
407  }
408  return jSmall;
409 }
410 
412 {
413  return new BandMatrix(*this);
414 }
415 
416 doublereal* BandMatrix::ptrColumn(size_t j)
417 {
418  return m_colPtrs[j];
419 }
420 
421 doublereal* const* BandMatrix::colPts()
422 {
423  return &(m_colPtrs[0]);
424 }
425 
427 {
428  warn_deprecated("BandMatrix::copyData", "To be removed after Cantera 2.2.");
429  m_factored = false;
430  size_t n = sizeof(doublereal) * m_n * (2 *m_kl + m_ku + 1);
431  GeneralMatrix* yyPtr = const_cast<GeneralMatrix*>(&y);
432  (void) memcpy(DATA_PTR(data), yyPtr->ptrColumn(0), n);
433 }
434 
435 void BandMatrix::useFactorAlgorithm(int fAlgorithm)
436 {
437 // useQR_ = fAlgorithm;
438 }
439 
440 
441 }
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)=0
Return a pointer to the top of column j, columns are assumed to be contiguous in memory.
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:39
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
Various templated functions that carry out common vector operations (see Templated Utility Functions)...
size_t nSuperDiagonals() const
Number of superdiagonals.
Definition: BandMatrix.cpp:183
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:165
int m_factored
Indicates whether the matrix is factored.
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
void warn_deprecated(const std::string &method, const std::string &extra)
Print a warning indicating that method is deprecated.
Definition: global.cpp:78
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
void writelogf(const char *fmt,...)
Write a formatted message to the screen.
Definition: global.cpp:38
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
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:99
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
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
Contains declarations for string manipulation functions within Cantera.
int factor()
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
Definition: BandMatrix.cpp:233
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
Definition: ct_defs.h:36
virtual void mult(const doublereal *b, doublereal *prod) const
Multiply A*b and write result to prod.
Definition: BandMatrix.cpp:198
virtual void zero()
Zero the matrix elements.
Definition: BandMatrix.cpp:113
std::vector< doublereal * > m_colPtrs
Vector of column pointers.
Definition: BandMatrix.h:336
GeneralMatrix & operator=(const GeneralMatrix &right)
Assignment operator.
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
Declarations for the class BandMatrix which is a child class of GeneralMatrix for banded matrices han...
A class for banded matrices, involving matrix inversion processes.
Definition: BandMatrix.h:33