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