Cantera  2.5.1
BandMatrix.cpp
Go to the documentation of this file.
1 //! @file BandMatrix.cpp Banded matrices.
2 
3 // This file is part of Cantera. See License.txt in the top-level directory or
4 // at https://cantera.org/license.txt for license and copyright information.
5 
9 
10 #if CT_USE_LAPACK
12 #else
13  #if CT_SUNDIALS_USE_LAPACK
14  #if CT_SUNDIALS_VERSION >= 30
15  #include "sunlinsol/sunlinsol_lapackband.h"
16  #else
17  #include "cvodes/cvodes_lapack.h"
18  #endif
19  #else
20  #if CT_SUNDIALS_VERSION >= 30
21  #include "sunlinsol/sunlinsol_band.h"
22  #else
23  #include "cvodes/cvodes_dense.h"
24  #include "cvodes/cvodes_band.h"
25  #endif
26  #endif
27 #endif
28 
29 #include <cstring>
30 #include <fstream>
31 
32 using namespace std;
33 
34 namespace Cantera
35 {
36 
37 // pImpl wrapper class for vector of Sundials index types to avoid needing to
38 // include Sundials headers in BandMatrix.h
39 struct BandMatrix::PivData {
40 #if CT_USE_LAPACK
41  vector_int data;
42 #elif CT_SUNDIALS_VERSION >= 30
43  std::vector<sunindextype> data;
44 #elif CT_SUNDIALS_VERSION >= 25
45  std::vector<long int> data;
46 #else
47  std::vector<int> data;
48 #endif
49 };
50 
51 BandMatrix::BandMatrix() :
52  m_n(0),
53  m_kl(0),
54  m_ku(0),
55  m_zero(0.0),
56  m_ipiv{new PivData()},
57  m_info(0)
58 {
59 }
60 
61 BandMatrix::~BandMatrix()
62 {
63  // Needs to be defined here so m_ipiv can be deleted
64 }
65 
66 BandMatrix::BandMatrix(size_t n, size_t kl, size_t ku, doublereal v) :
67  m_n(n),
68  m_kl(kl),
69  m_ku(ku),
70  m_zero(0.0),
71  m_ipiv{new PivData()},
72  m_info(0)
73 {
74  data.resize(n*(2*kl + ku + 1));
75  ludata.resize(n*(2*kl + ku + 1));
76  fill(data.begin(), data.end(), v);
77  fill(ludata.begin(), ludata.end(), 0.0);
78  m_ipiv->data.resize(m_n);
79  m_colPtrs.resize(n);
80  m_lu_col_ptrs.resize(n);
81  size_t ldab = (2*kl + ku + 1);
82  for (size_t j = 0; j < n; j++) {
83  m_colPtrs[j] = &data[ldab * j];
84  m_lu_col_ptrs[j] = &ludata[ldab * j];
85  }
86 }
87 
89  GeneralMatrix(y),
90  m_n(0),
91  m_kl(0),
92  m_ku(0),
93  m_zero(0.0),
94  m_ipiv{new PivData()},
95  m_info(y.m_info)
96 {
97  m_n = y.m_n;
98  m_kl = y.m_kl;
99  m_ku = y.m_ku;
100  data = y.data;
101  ludata = y.ludata;
102  m_ipiv->data = y.m_ipiv->data;
103  m_colPtrs.resize(m_n);
104  m_lu_col_ptrs.resize(m_n);
105  size_t ldab = (2 *m_kl + m_ku + 1);
106  for (size_t j = 0; j < m_n; j++) {
107  m_colPtrs[j] = &data[ldab * j];
108  m_lu_col_ptrs[j] = &ludata[ldab * j];
109  }
110 }
111 
112 BandMatrix& BandMatrix::operator=(const BandMatrix& y)
113 {
114  if (&y == this) {
115  return *this;
116  }
117  GeneralMatrix::operator=(y);
118  m_n = y.m_n;
119  m_kl = y.m_kl;
120  m_ku = y.m_ku;
121  m_ipiv->data = y.m_ipiv->data;
122  data = y.data;
123  ludata = y.ludata;
124  m_colPtrs.resize(m_n);
125  m_lu_col_ptrs.resize(m_n);
126  size_t ldab = (2 * m_kl + m_ku + 1);
127  for (size_t j = 0; j < m_n; j++) {
128  m_colPtrs[j] = &data[ldab * j];
129  m_lu_col_ptrs[j] = &ludata[ldab * j];
130  }
131  m_info = y.m_info;
132  return *this;
133 }
134 
135 void BandMatrix::resize(size_t n, size_t kl, size_t ku, doublereal v)
136 {
137  m_n = n;
138  m_kl = kl;
139  m_ku = ku;
140  data.resize(n*(2*kl + ku + 1));
141  ludata.resize(n*(2*kl + ku + 1));
142  m_ipiv->data.resize(m_n);
143  fill(data.begin(), data.end(), v);
144  m_colPtrs.resize(m_n);
145  m_lu_col_ptrs.resize(m_n);
146  size_t ldab = (2 * m_kl + m_ku + 1);
147  for (size_t j = 0; j < n; j++) {
148  m_colPtrs[j] = &data[ldab * j];
149  m_lu_col_ptrs[j] = &ludata[ldab * j];
150  }
151  m_factored = false;
152 }
153 
154 void BandMatrix::bfill(doublereal v)
155 {
156  std::fill(data.begin(), data.end(), v);
157  m_factored = false;
158 }
159 
161 {
162  std::fill(data.begin(), data.end(), 0.0);
163  m_factored = false;
164 }
165 
166 doublereal& BandMatrix::operator()(size_t i, size_t j)
167 {
168  return value(i,j);
169 }
170 
171 doublereal BandMatrix::operator()(size_t i, size_t j) const
172 {
173  return value(i,j);
174 }
175 
176 doublereal& BandMatrix::value(size_t i, size_t j)
177 {
178  m_factored = false;
179  if (i + m_ku < j || i > j + m_kl) {
180  return m_zero;
181  }
182  return data[index(i,j)];
183 }
184 
185 doublereal BandMatrix::value(size_t i, size_t j) const
186 {
187  if (i + m_ku < j || i > j + m_kl) {
188  return 0.0;
189  }
190  return data[index(i,j)];
191 }
192 
193 size_t BandMatrix::index(size_t i, size_t j) const
194 {
195  return (2*m_kl + m_ku)*j + m_kl + m_ku + i;
196 }
197 
198 doublereal BandMatrix::_value(size_t i, size_t j) const
199 {
200  return data[index(i,j)];
201 }
202 
203 size_t BandMatrix::nRows() const
204 {
205  return m_n;
206 }
207 
208 size_t BandMatrix::nColumns() const
209 {
210  return m_n;
211 }
212 
214 {
215  return m_kl;
216 }
217 
219 {
220  return m_ku;
221 }
222 
223 size_t BandMatrix::ldim() const
224 {
225  return 2*m_kl + m_ku + 1;
226 }
227 
228 void BandMatrix::mult(const doublereal* b, doublereal* prod) const
229 {
230  for (size_t m = 0; m < m_n; m++) {
231  double sum = 0.0;
232  size_t start = (m >= m_kl) ? m - m_kl : 0;
233  size_t stop = std::min(m + m_ku + 1, m_n);
234  for (size_t j = start; j < stop; j++) {
235  sum += _value(m,j) * b[j];
236  }
237  prod[m] = sum;
238  }
239 }
240 
241 void BandMatrix::leftMult(const doublereal* const b, doublereal* const prod) const
242 {
243  for (size_t n = 0; n < m_n; n++) {
244  double sum = 0.0;
245  size_t start = (n >= m_ku) ? n - m_ku : 0;
246  size_t stop = std::min(n + m_kl + 1, m_n);
247  for (size_t i = start; i < stop; i++) {
248  sum += _value(i,n) * b[i];
249  }
250  prod[n] = sum;
251  }
252 }
253 
255 {
256  ludata = data;
257 #if CT_USE_LAPACK
258  ct_dgbtrf(nRows(), nColumns(), nSubDiagonals(), nSuperDiagonals(),
259  ludata.data(), ldim(), m_ipiv->data.data(), m_info);
260 #else
261  long int nu = static_cast<long int>(nSuperDiagonals());
262  long int nl = static_cast<long int>(nSubDiagonals());
263  long int smu = nu + nl;
264  m_info = bandGBTRF(m_lu_col_ptrs.data(), static_cast<long int>(nColumns()),
265  nu, nl, smu, m_ipiv->data.data());
266 #endif
267  if (m_info != 0) {
268  throw Cantera::CanteraError("BandMatrix::factor",
269  "Factorization failed with DGBTRF error code {}.", m_info);
270  }
271  m_factored = true;
272  return m_info;
273 }
274 
275 int BandMatrix::solve(const doublereal* const b, doublereal* const x)
276 {
277  copy(b, b + m_n, x);
278  return solve(x);
279 }
280 
281 int BandMatrix::solve(doublereal* b, size_t nrhs, size_t ldb)
282 {
283  if (!m_factored) {
284  factor();
285  }
286  if (ldb == 0) {
287  ldb = nColumns();
288  }
289 #if CT_USE_LAPACK
290  ct_dgbtrs(ctlapack::NoTranspose, nColumns(), nSubDiagonals(),
291  nSuperDiagonals(), nrhs, ludata.data(), ldim(),
292  m_ipiv->data.data(), b, ldb, m_info);
293 #else
294  long int nu = static_cast<long int>(nSuperDiagonals());
295  long int nl = static_cast<long int>(nSubDiagonals());
296  long int smu = nu + nl;
297  double** a = m_lu_col_ptrs.data();
298  bandGBTRS(a, static_cast<long int>(nColumns()), smu, nl, m_ipiv->data.data(), b);
299  m_info = 0;
300 #endif
301 
302  if (m_info != 0) {
303  throw Cantera::CanteraError("BandMatrix::solve",
304  "Linear solve failed with DGBTRS error code {}.", m_info);
305  }
306  return m_info;
307 }
308 
309 vector_fp::iterator BandMatrix::begin()
310 {
311  m_factored = false;
312  return data.begin();
313 }
314 
315 vector_fp::iterator BandMatrix::end()
316 {
317  m_factored = false;
318  return data.end();
319 }
320 
321 vector_fp::const_iterator BandMatrix::begin() const
322 {
323  return data.begin();
324 }
325 
326 vector_fp::const_iterator BandMatrix::end() const
327 {
328  return data.end();
329 }
330 
331 ostream& operator<<(ostream& s, const BandMatrix& m)
332 {
333  for (size_t i = 0; i < m.nRows(); i++) {
334  s << m(i, 0);
335  for (size_t j = 1; j < m.nColumns(); j++) {
336  s << ", " << m(i,j);
337  }
338  s << endl;
339  }
340  return s;
341 }
342 
343 doublereal BandMatrix::rcond(doublereal a1norm)
344 {
345  iwork_.resize(m_n);
346  work_.resize(3 * m_n);
347 
348  if (m_factored != 1) {
349  throw CanteraError("BandMatrix::rcond", "matrix isn't factored correctly");
350  }
351 
352 #if CT_USE_LAPACK
353  size_t ldab = (2 *m_kl + m_ku + 1);
354  int rinfo = 0;
355  double rcond = ct_dgbcon('1', m_n, m_kl, m_ku, ludata.data(),
356  ldab, m_ipiv->data.data(), a1norm, work_.data(), iwork_.data(), rinfo);
357  if (rinfo != 0) {
358  throw CanteraError("BandMatrix::rcond", "DGBCON returned INFO = {}", rinfo);
359  }
360  return rcond;
361 #else
362  throw CanteraError("BandMatrix::rcond", "not implemented when LAPACK is missing");
363 #endif
364 }
365 
367 {
368  return 0;
369 }
370 
371 doublereal BandMatrix::oneNorm() const
372 {
373  double value = 0.0;
374  for (size_t j = 0; j < m_n; j++) {
375  double sum = 0.0;
376  size_t start = (j >= m_ku) ? j - m_ku : 0;
377  size_t stop = std::min(j + m_kl + 1, m_n);
378  for (size_t i = start; i < stop; i++) {
379  sum += std::abs(_value(i,j));
380  }
381  value = std::max(sum, value);
382  }
383  return value;
384 }
385 
386 size_t BandMatrix::checkRows(doublereal& valueSmall) const
387 {
388  valueSmall = 1.0E300;
389  size_t iSmall = npos;
390  for (size_t i = 0; i < m_n; i++) {
391  double valueS = 0.0;
392  size_t start = (i > m_kl) ? i - m_kl : 0;
393  size_t stop = std::min(i + m_ku + 1, m_n);
394  for (size_t j = start; j < stop; j++) {
395  valueS = std::max(fabs(_value(i,j)), valueS);
396  }
397  if (valueS < valueSmall) {
398  iSmall = i;
399  valueSmall = valueS;
400  if (valueSmall == 0.0) {
401  return iSmall;
402  }
403  }
404  }
405  return iSmall;
406 }
407 
408 size_t BandMatrix::checkColumns(doublereal& valueSmall) const
409 {
410  valueSmall = 1.0E300;
411  size_t jSmall = npos;
412  for (size_t j = 0; j < m_n; j++) {
413  double valueS = 0.0;
414  size_t start = (j > m_ku) ? j - m_ku : 0;
415  size_t stop = std::min(j + m_kl + 1, m_n);
416  for (size_t i = start; i < stop; i++) {
417  valueS = std::max(fabs(_value(i,j)), valueS);
418  }
419  if (valueS < valueSmall) {
420  jSmall = j;
421  valueSmall = valueS;
422  if (valueSmall == 0.0) {
423  return jSmall;
424  }
425  }
426  }
427  return jSmall;
428 }
429 
430 doublereal* BandMatrix::ptrColumn(size_t j)
431 {
432  return m_colPtrs[j];
433 }
434 
435 doublereal* const* BandMatrix::colPts()
436 {
437  return &m_colPtrs[0];
438 }
439 
440 }
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:35
virtual void zero()
Zero the matrix elements.
Definition: BandMatrix.cpp:160
virtual doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j.
Definition: BandMatrix.cpp:430
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:198
size_t m_ku
Number of super diagonals of the matrix.
Definition: BandMatrix.h:277
size_t m_kl
Number of subdiagonals of the matrix.
Definition: BandMatrix.h:274
vector_fp ludata
Factorized data.
Definition: BandMatrix.h:268
virtual size_t checkRows(doublereal &valueSmall) const
Check to see if we have any zero rows in the Jacobian.
Definition: BandMatrix.cpp:386
std::vector< doublereal * > m_colPtrs
Vector of column pointers.
Definition: BandMatrix.h:288
size_t nSubDiagonals() const
Number of subdiagonals.
Definition: BandMatrix.cpp:213
size_t ldim() const
Return the number of rows of storage needed for the band storage.
Definition: BandMatrix.cpp:223
int factor()
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
Definition: BandMatrix.cpp:254
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:193
virtual void mult(const doublereal *b, doublereal *prod) const
Multiply A*b and write result to prod.
Definition: BandMatrix.cpp:228
virtual size_t nRows() const
Return the number of rows in the matrix.
Definition: BandMatrix.cpp:203
virtual doublereal rcond(doublereal a1norm)
Returns an estimate of the inverse of the condition number for the matrix.
Definition: BandMatrix.cpp:343
vector_int iwork_
Extra work array needed - size = n.
Definition: BandMatrix.h:292
virtual vector_fp::iterator begin()
Returns an iterator for the start of the band storage data.
Definition: BandMatrix.cpp:309
std::unique_ptr< PivData > m_ipiv
Pivot vector.
Definition: BandMatrix.h:282
void bfill(doublereal v=0.0)
Fill or zero the matrix.
Definition: BandMatrix.cpp:154
BandMatrix()
Base Constructor.
Definition: BandMatrix.cpp:51
size_t nColumns() const
Number of columns.
Definition: BandMatrix.cpp:208
doublereal m_zero
value of zero
Definition: BandMatrix.h:280
doublereal & operator()(size_t i, size_t j)
Index into the (i,j) element.
Definition: BandMatrix.cpp:166
size_t nSuperDiagonals() const
Number of superdiagonals.
Definition: BandMatrix.cpp:218
vector_fp work_
Extra dp work array needed - size = 3n.
Definition: BandMatrix.h:295
int solve(const doublereal *const b, doublereal *const x)
Solve the matrix problem Ax = b.
Definition: BandMatrix.cpp:275
virtual int factorAlgorithm() const
Returns the factor algorithm used.
Definition: BandMatrix.cpp:366
vector_fp::iterator end()
Returns an iterator for the end of the band storage data.
Definition: BandMatrix.cpp:315
void resize(size_t n, size_t kl, size_t ku, doublereal v=0.0)
Resize the matrix problem.
Definition: BandMatrix.cpp:135
virtual size_t checkColumns(doublereal &valueSmall) const
Check to see if we have any zero columns in the Jacobian.
Definition: BandMatrix.cpp:408
virtual void leftMult(const doublereal *const b, doublereal *const prod) const
Multiply b*A and write result to prod.
Definition: BandMatrix.cpp:241
virtual doublereal *const * colPts()
Return a vector of const pointers to the columns.
Definition: BandMatrix.cpp:435
doublereal & value(size_t i, size_t j)
Return a changeable reference to element (i,j).
Definition: BandMatrix.cpp:176
size_t m_n
Number of rows and columns of the matrix.
Definition: BandMatrix.h:271
vector_fp data
Matrix data.
Definition: BandMatrix.h:261
virtual doublereal oneNorm() const
Returns the one norm of the matrix.
Definition: BandMatrix.cpp:371
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
Generic matrix.
Definition: GeneralMatrix.h:22
int m_factored
Indicates whether the matrix is factored.
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:188
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:182
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector operations (see Templated Utility Functions)...