Cantera  2.1.2
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"
14 #include "cantera/base/global.h"
15 
16 #include <cstring>
17 #include <fstream>
18 
19 using namespace std;
20 
21 namespace Cantera
22 {
23 
24 BandMatrix::BandMatrix() :
25  GeneralMatrix(1),
26  m_factored(false),
27  m_n(0),
28  m_kl(0),
29  m_ku(0),
30  m_zero(0.0)
31 {
32  data.clear();
33  ludata.clear();
34 }
35 
36 BandMatrix::BandMatrix(size_t n, size_t kl, size_t ku, doublereal v) :
37  GeneralMatrix(1),
38  m_factored(false),
39  m_n(n),
40  m_kl(kl),
41  m_ku(ku),
42  m_zero(0.0)
43 {
44  data.resize(n*(2*kl + ku + 1));
45  ludata.resize(n*(2*kl + ku + 1));
46  fill(data.begin(), data.end(), v);
47  fill(ludata.begin(), ludata.end(), 0.0);
48  m_ipiv.resize(m_n);
49  m_colPtrs.resize(n);
50  size_t ldab = (2*kl + ku + 1);
51  for (size_t j = 0; j < n; j++) {
52  m_colPtrs[j] = &(data[ldab * j]);
53  }
54 }
55 
57  GeneralMatrix(1),
58  m_factored(false),
59  m_n(0),
60  m_kl(0),
61  m_ku(0),
62  m_zero(0.0)
63 {
64  m_n = y.m_n;
65  m_kl = y.m_kl;
66  m_ku = y.m_ku;
67  data = y.data;
68  ludata = y.ludata;
70  m_ipiv = y.m_ipiv;
71  m_colPtrs.resize(m_n);
72  size_t ldab = (2 *m_kl + m_ku + 1);
73  for (size_t j = 0; j < m_n; j++) {
74  m_colPtrs[j] = &(data[ldab * j]);
75  }
76 }
77 
79 {
80  if (&y == this) {
81  return *this;
82  }
84  m_n = y.m_n;
85  m_kl = y.m_kl;
86  m_ku = y.m_ku;
87  m_ipiv = y.m_ipiv;
88  data = y.data;
89  ludata = y.ludata;
91  m_colPtrs.resize(m_n);
92  size_t ldab = (2 * m_kl + m_ku + 1);
93  for (size_t j = 0; j < m_n; j++) {
94  m_colPtrs[j] = &(data[ldab * j]);
95  }
96  return *this;
97 }
98 
99 void BandMatrix::resize(size_t n, size_t kl, size_t ku, doublereal v)
100 {
101  m_n = n;
102  m_kl = kl;
103  m_ku = ku;
104  data.resize(n*(2*kl + ku + 1));
105  ludata.resize(n*(2*kl + ku + 1));
106  m_ipiv.resize(m_n);
107  fill(data.begin(), data.end(), v);
108  m_colPtrs.resize(m_n);
109  size_t ldab = (2 * m_kl + m_ku + 1);
110  for (size_t j = 0; j < n; j++) {
111  m_colPtrs[j] = &(data[ldab * j]);
112  }
113  m_factored = false;
114 }
115 
116 void BandMatrix::bfill(doublereal v)
117 {
118  std::fill(data.begin(), data.end(), v);
119  m_factored = false;
120 }
121 
123 {
124  std::fill(data.begin(), data.end(), 0.0);
125  m_factored = false;
126 }
127 
128 doublereal& BandMatrix::operator()(size_t i, size_t j)
129 {
130  return value(i,j);
131 }
132 
133 doublereal BandMatrix::operator()(size_t i, size_t j) const
134 {
135  return value(i,j);
136 }
137 
138 doublereal& BandMatrix::value(size_t i, size_t j)
139 {
140  m_factored = false;
141  if (i + m_ku < j || i > j + m_kl) {
142  return m_zero;
143  }
144  return data[index(i,j)];
145 }
146 
147 doublereal BandMatrix::value(size_t i, size_t j) const
148 {
149  if (i + m_ku < j || i > j + m_kl) {
150  return 0.0;
151  }
152  return data[index(i,j)];
153 }
154 
155 size_t BandMatrix::index(size_t i, size_t j) const
156 {
157  int jj = static_cast<int>(j);
158  int ii = static_cast<int>(i);
159  size_t rw = (int) m_kl + (int) m_ku + (int) ii - jj;
160  return (2*m_kl + m_ku + 1)*j + rw;
161 }
162 
163 doublereal BandMatrix::_value(size_t i, size_t j) const
164 {
165  return data[index(i,j)];
166 }
167 
168 size_t BandMatrix::nRows() const
169 {
170  return m_n;
171 }
172 
173 size_t BandMatrix::nRowsAndStruct(size_t* const iStruct) const
174 {
175  if (iStruct) {
176  iStruct[0] = m_kl;
177  iStruct[1] = m_ku;
178  }
179  return m_n;
180 }
181 //====================================================================================================================
182 // Number of columns
183 size_t BandMatrix::nColumns() const
184 {
185  return m_n;
186 }
187 //====================================================================================================================
188 // Number of subdiagonals
190 {
191  return m_kl;
192 }
193 //====================================================================================================================
194 // Number of superdiagonals
196 {
197  return m_ku;
198 }
199 //====================================================================================================================
200 size_t BandMatrix::ldim() const
201 {
202  return 2*m_kl + m_ku + 1;
203 }
204 //====================================================================================================================
206 {
207  return m_ipiv;
208 }
209 //====================================================================================================================
210 /*
211  * Multiply A*b and write result to \c prod.
212  */
213 void BandMatrix::mult(const doublereal* b, doublereal* prod) const
214 {
215  int kl = static_cast<int>(m_kl);
216  int ku = static_cast<int>(m_ku);
217  int nr = static_cast<int>(nRows());
218  doublereal sum = 0.0;
219  for (int m = 0; m < nr; m++) {
220  sum = 0.0;
221  for (int j = m - kl; j <= m + ku; j++) {
222  if (j >= 0 && j < (int) m_n) {
223  sum += _value(m,j) * b[j];
224  }
225  }
226  prod[m] = sum;
227  }
228 }
229 
230 void BandMatrix::leftMult(const doublereal* const b, doublereal* const prod) const
231 {
232  int kl = static_cast<int>(m_kl);
233  int ku = static_cast<int>(m_ku);
234  int nc = static_cast<int>(nColumns());
235  doublereal sum = 0.0;
236  for (int n = 0; n < nc; n++) {
237  sum = 0.0;
238  for (int i = n - ku; i <= n + kl; i++) {
239  if (i >= 0 && i < (int) m_n) {
240  size_t ii = i;
241  sum += _value(ii,n) * b[ii];
242  }
243  }
244  prod[n] = sum;
245  }
246 }
247 
249 {
250  int info=0;
251  copy(data.begin(), data.end(), ludata.begin());
252  ct_dgbtrf(nRows(), nColumns(), nSubDiagonals(), nSuperDiagonals(),
253  DATA_PTR(ludata), ldim(), DATA_PTR(ipiv()), info);
254 
255  // if info = 0, LU decomp succeeded.
256  if (info == 0) {
257  m_factored = true;
258  } else {
259  m_factored = false;
260  ofstream fout("bandmatrix.csv");
261  fout << *this << endl;
262  fout.close();
263  }
264  return info;
265 }
266 
267 int BandMatrix::solve(const doublereal* const b, doublereal* const x)
268 {
269  copy(b, b + m_n, x);
270  return solve(x);
271 }
272 
273 int BandMatrix::solve(doublereal* b)
274 {
275  int info = 0;
276  if (!m_factored) {
277  info = factor();
278  }
279  if (info == 0)
280  ct_dgbtrs(ctlapack::NoTranspose, nColumns(), nSubDiagonals(),
282  DATA_PTR(ipiv()), b, nColumns(), info);
283 
284  // error handling
285  if (info != 0) {
286  ofstream fout("bandmatrix.csv");
287  fout << *this << endl;
288  fout.close();
289  }
290  return info;
291 }
292 
293 vector_fp::iterator BandMatrix::begin()
294 {
295  m_factored = false;
296  return data.begin();
297 }
298 
299 vector_fp::iterator BandMatrix::end()
300 {
301  m_factored = false;
302  return data.end();
303 }
304 
305 vector_fp::const_iterator BandMatrix::begin() const
306 {
307  return data.begin();
308 }
309 
310 vector_fp::const_iterator BandMatrix::end() const
311 {
312  return data.end();
313 }
314 
315 ostream& operator<<(ostream& s, const BandMatrix& m)
316 {
317  size_t nr = m.nRows();
318  size_t nc = m.nColumns();
319  for (size_t i = 0; i < nr; i++) {
320  for (size_t j = 0; j < nc; j++) {
321  s << m(i,j) << ", ";
322  }
323  s << endl;
324  }
325  return s;
326 }
327 
328 void BandMatrix::err(const std::string& msg) const
329 {
330  throw CanteraError("BandMatrix() unimplemented function", msg);
331 }
332 
334 {
335  factor();
336  return 0;
337 }
338 
340 {
341  double a1norm = oneNorm();
342  return rcond(a1norm);
343 }
344 
345 doublereal BandMatrix::rcond(doublereal a1norm)
346 {
347  int printLevel = 0;
348  int useReturnErrorCode = 0;
349  if (iwork_.size() < m_n) {
350  iwork_.resize(m_n);
351  }
352  if (work_.size() < 3 * m_n) {
353  work_.resize(3 * m_n);
354  }
355  doublereal rcond = 0.0;
356  if (m_factored != 1) {
357  throw CanteraError("BandMatrix::rcond()", "matrix isn't factored correctly");
358  }
359 
360  // doublereal anorm = oneNorm();
361  size_t ldab = (2 *m_kl + m_ku + 1);
362  int rinfo = 0;
363  rcond = ct_dgbcon('1', m_n, m_kl, m_ku, DATA_PTR(ludata), ldab, DATA_PTR(m_ipiv), a1norm, DATA_PTR(work_),
364  DATA_PTR(iwork_), rinfo);
365  if (rinfo != 0) {
366  if (printLevel) {
367  writelogf("BandMatrix::rcond(): DGBCON returned INFO = %d\n", rinfo);
368  }
369  if (! useReturnErrorCode) {
370  throw CanteraError("BandMatrix::rcond()", "DGBCON returned INFO = " + int2str(rinfo));
371  }
372  }
373  return rcond;
374 }
375 
376 void BandMatrix::useFactorAlgorithm(int fAlgorithm)
377 {
378  // QR algorithm isn't implemented for banded matrix.
379 }
380 
382 {
383  return 0;
384 }
385 
386 doublereal BandMatrix::oneNorm() const
387 {
388  int ku = static_cast<int>(m_ku);
389  int kl = static_cast<int>(m_kl);
390  doublereal value = 0.0;
391  for (int j = 0; j < (int) m_n; j++) {
392  doublereal sum = 0.0;
393  doublereal* colP = m_colPtrs[j];
394  for (int i = j - ku; i <= j + kl; i++) {
395  sum += fabs(colP[kl + ku + i - j]);
396  }
397  if (sum > value) {
398  value = sum;
399  }
400  }
401  return value;
402 }
403 
404 size_t BandMatrix::checkRows(doublereal& valueSmall) const
405 {
406  valueSmall = 1.0E300;
407  size_t iSmall = npos;
408  double vv;
409  for (int i = 0; i < (int) m_n; i++) {
410  double valueS = 0.0;
411  for (int j = i - (int) m_kl; j <= i + (int) m_ku; j++) {
412  if (j >= 0 && j < (int) m_n) {
413  vv = fabs(value(i,j));
414  if (vv > valueS) {
415  valueS = vv;
416  }
417  }
418  }
419  if (valueS < valueSmall) {
420  iSmall = i;
421  valueSmall = valueS;
422  if (valueSmall == 0.0) {
423  return iSmall;
424  }
425  }
426  }
427  return iSmall;
428 }
429 
430 size_t BandMatrix::checkColumns(doublereal& valueSmall) const
431 {
432  valueSmall = 1.0E300;
433  size_t jSmall = npos;
434  double vv;
435  for (int j = 0; j < (int) m_n; j++) {
436  double valueS = 0.0;
437  for (int i = j - (int) m_ku; i <= j + (int) m_kl; i++) {
438  if (i >= 0 && i < (int) m_n) {
439  vv = fabs(value(i,j));
440  if (vv > valueS) {
441  valueS = vv;
442  }
443  }
444  }
445  if (valueS < valueSmall) {
446  jSmall = j;
447  valueSmall = valueS;
448  if (valueSmall == 0.0) {
449  return jSmall;
450  }
451  }
452  }
453  return jSmall;
454 }
455 
457 {
458  return new BandMatrix(*this);
459 }
460 
462 {
463  return m_factored;
464 }
465 
466 doublereal* BandMatrix::ptrColumn(size_t j)
467 {
468  return m_colPtrs[j];
469 }
470 
471 doublereal* const* BandMatrix::colPts()
472 {
473  return &(m_colPtrs[0]);
474 }
475 
477 {
478  m_factored = false;
479  size_t n = sizeof(doublereal) * m_n * (2 *m_kl + m_ku + 1);
480  GeneralMatrix* yyPtr = const_cast<GeneralMatrix*>(&y);
481  (void) memcpy(DATA_PTR(data), yyPtr->ptrColumn(0), n);
482 }
483 
485 {
486  m_factored = 0;
487 }
488 
489 }
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)=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:40
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
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:173
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
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
This file contains definitions for utility functions and text for modules, inputfiles, logs, textlogs, HTML_logs (see Input File Handling, Diagnostic Output, Writing messages to the screen and Writing HTML Logfiles).
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 writelogf(const char *fmt,...)
Write a formatted message to the screen.
Definition: global.cpp:48
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
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:68
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
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
Contains declarations for string manipulation functions within Cantera.
int factor()
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
Definition: BandMatrix.cpp:248
#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:213
virtual void zero()
Zero the matrix elements.
Definition: BandMatrix.cpp:122
std::vector< doublereal * > m_colPtrs
Vector of column pointers.
Definition: BandMatrix.h:365
GeneralMatrix & operator=(const GeneralMatrix &right)
Assignment operator.
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
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:37