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