Cantera  4.0.0a1
Loading...
Searching...
No Matches
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 #include "sunlinsol/sunlinsol_band.h"
14#endif
15
16#include <cstring>
17#include <fstream>
18
19using namespace std;
20
21namespace Cantera
22{
23
24// pImpl wrapper class for vector of Sundials index types to avoid needing to
25// include Sundials headers in BandMatrix.h
27#if CT_USE_LAPACK
28 vector<int> data;
29#else
30 vector<sunindextype> data;
31#endif
32};
33
35 m_ipiv{new PivData()}
36{
37}
38
39BandMatrix::~BandMatrix()
40{
41 // Needs to be defined here so m_ipiv can be deleted
42}
43
44BandMatrix::BandMatrix(size_t n, size_t kl, size_t ku, double v) :
45 m_n(n),
46 m_kl(kl),
47 m_ku(ku),
48 m_ipiv{new PivData()}
49{
50 data.resize(n*(2*kl + ku + 1));
51 ludata.resize(n*(2*kl + ku + 1));
52 fill(data.begin(), data.end(), v);
53 fill(ludata.begin(), ludata.end(), 0.0);
54 m_ipiv->data.resize(m_n);
55 m_lu_col_ptrs.resize(n);
56 size_t ldab = (2*kl + ku + 1);
57 for (size_t j = 0; j < n; j++) {
58 m_lu_col_ptrs[j] = &ludata[ldab * j];
59 }
60}
61
64 data(y.data),
65 ludata(y.ludata),
66 m_n(y.m_n),
67 m_kl(y.m_kl),
68 m_ku(y.m_ku),
69 m_ipiv{new PivData()},
70 m_info(y.m_info)
71{
72 m_ipiv->data = y.m_ipiv->data;
73 m_lu_col_ptrs.resize(m_n);
74 size_t ldab = (2*m_kl + m_ku + 1);
75 for (size_t j = 0; j < m_n; j++) {
76 m_lu_col_ptrs[j] = &ludata[ldab * j];
77 }
78}
79
80BandMatrix& BandMatrix::operator=(const BandMatrix& y)
81{
82 if (&y == this) {
83 return *this;
84 }
85 GeneralMatrix::operator=(y);
86 m_n = y.m_n;
87 m_kl = y.m_kl;
88 m_ku = y.m_ku;
89 m_ipiv->data = y.m_ipiv->data;
90 data = y.data;
91 ludata = y.ludata;
92 m_lu_col_ptrs.resize(m_n);
93 size_t ldab = (2*m_kl + m_ku + 1);
94 for (size_t j = 0; j < m_n; j++) {
95 m_lu_col_ptrs[j] = &ludata[ldab * j];
96 }
97 m_info = y.m_info;
98 return *this;
99}
100
101void BandMatrix::resize(size_t n, size_t kl, size_t ku, double v)
102{
103 m_n = n;
104 m_kl = kl;
105 m_ku = ku;
106 data.resize(n*(2*kl + ku + 1));
107 ludata.resize(n*(2*kl + ku + 1));
108 m_lu_col_ptrs.resize(m_n);
109 size_t ldab = (2*m_kl + m_ku + 1);
110 for (size_t j = 0; j < m_n; j++) {
111 m_lu_col_ptrs[j] = &ludata[ldab * j];
112 }
113 m_ipiv->data.resize(m_n);
114 fill(data.begin(), data.end(), v);
115 m_factored = false;
116}
117
118void BandMatrix::bfill(double v)
119{
120 std::fill(data.begin(), data.end(), v);
121 m_factored = false;
122}
123
125{
126 std::fill(data.begin(), data.end(), 0.0);
127 m_factored = false;
128}
129
130double& BandMatrix::operator()(size_t i, size_t j)
131{
132 return value(i,j);
133}
134
135double BandMatrix::operator()(size_t i, size_t j) const
136{
137 return value(i,j);
138}
139
140double& BandMatrix::value(size_t i, size_t j)
141{
142 m_factored = false;
143 if (i + m_ku < j || i > j + m_kl) {
144 return m_zero;
145 }
146 return data[index(i,j)];
147}
148
149double BandMatrix::value(size_t i, size_t j) const
150{
151 if (i + m_ku < j || i > j + m_kl) {
152 return 0.0;
153 }
154 return data[index(i,j)];
155}
156
157size_t BandMatrix::index(size_t i, size_t j) const
158{
159 return (2*m_kl + m_ku)*j + m_kl + m_ku + i;
160}
161
162double BandMatrix::_value(size_t i, size_t j) const
163{
164 return data[index(i,j)];
165}
166
167size_t BandMatrix::nRows() const
168{
169 return m_n;
170}
171
173{
174 return m_n;
175}
176
178{
179 return m_kl;
180}
181
183{
184 return m_ku;
185}
186
187size_t BandMatrix::ldim() const
188{
189 return 2*m_kl + m_ku + 1;
190}
191
192void BandMatrix::mult(span<const double> b, span<double> prod) const
193{
194 checkArraySize("BandMatrix::mult", b.size(), m_n);
195 checkArraySize("BandMatrix::mult", prod.size(), m_n);
196 for (size_t m = 0; m < m_n; m++) {
197 double sum = 0.0;
198 size_t start = (m >= m_kl) ? m - m_kl : 0;
199 size_t stop = std::min(m + m_ku + 1, m_n);
200 for (size_t j = start; j < stop; j++) {
201 sum += _value(m,j) * b[j];
202 }
203 prod[m] = sum;
204 }
205}
206
207void BandMatrix::leftMult(span<const double> b, span<double> prod) const
208{
209 checkArraySize("BandMatrix::leftMult", b.size(), m_n);
210 checkArraySize("BandMatrix::leftMult", prod.size(), m_n);
211 for (size_t n = 0; n < m_n; n++) {
212 double sum = 0.0;
213 size_t start = (n >= m_ku) ? n - m_ku : 0;
214 size_t stop = std::min(n + m_kl + 1, m_n);
215 for (size_t i = start; i < stop; i++) {
216 sum += _value(i,n) * b[i];
217 }
218 prod[n] = sum;
219 }
220}
221
223{
224 ludata = data;
225#if CT_USE_LAPACK
226 ct_dgbtrf(nRows(), nColumns(), nSubDiagonals(), nSuperDiagonals(),
227 ludata.data(), ldim(), m_ipiv->data.data(), m_info);
228#else
229 long int nu = static_cast<long int>(nSuperDiagonals());
230 long int nl = static_cast<long int>(nSubDiagonals());
231 long int smu = nu + nl;
232 m_info = SUNDlsMat_bandGBTRF(m_lu_col_ptrs.data(),
233 static_cast<long int>(nColumns()),
234 nu, nl, smu, m_ipiv->data.data());
235#endif
236 if (m_info != 0) {
237 throw Cantera::CanteraError("BandMatrix::factor",
238 "Factorization failed with DGBTRF error code {}.", m_info);
239 }
240 m_factored = true;
241}
242
243void BandMatrix::solve(span<const double> b, span<double> x)
244{
245 checkArraySize("BandMatrix::solve", b.size(), m_n);
246 checkArraySize("BandMatrix::solve", x.size(), m_n);
247 copy(b.begin(), b.begin() + m_n, x.begin());
248 solve(x);
249}
250
251void BandMatrix::solve(span<double> b, size_t nrhs, size_t ldb)
252{
253 if (!m_factored) {
254 factor();
255 }
256 if (ldb == 0) {
257 ldb = nColumns();
258 }
259 checkArraySize("BandMatrix::solve", b.size(), nrhs * ldb);
260#if CT_USE_LAPACK
261 ct_dgbtrs(ctlapack::NoTranspose, nColumns(), nSubDiagonals(),
262 nSuperDiagonals(), nrhs, ludata.data(), ldim(),
263 m_ipiv->data.data(), b.data(), ldb, m_info);
264 if (m_info != 0) {
265 throw Cantera::CanteraError("BandMatrix::solve",
266 "Linear solve failed with DGBTRS error code {}.", m_info);
267 }
268#else
269 long int nu = static_cast<long int>(nSuperDiagonals());
270 long int nl = static_cast<long int>(nSubDiagonals());
271 long int smu = nu + nl;
272 double** a = m_lu_col_ptrs.data();
273 SUNDlsMat_bandGBTRS(a, static_cast<long int>(nColumns()), smu, nl,
274 m_ipiv->data.data(), b.data());
275 m_info = 0;
276#endif
277}
278
279ostream& operator<<(ostream& s, const BandMatrix& m)
280{
281 for (size_t i = 0; i < m.nRows(); i++) {
282 s << m(i, 0);
283 for (size_t j = 1; j < m.nColumns(); j++) {
284 s << ", " << m(i,j);
285 }
286 s << endl;
287 }
288 return s;
289}
290
291double BandMatrix::rcond(double a1norm)
292{
293 iwork_.resize(m_n);
294 work_.resize(3 * m_n);
295
296 if (m_factored != 1) {
297 throw CanteraError("BandMatrix::rcond", "matrix isn't factored correctly");
298 }
299
300#if CT_USE_LAPACK
301 size_t ldab = (2 *m_kl + m_ku + 1);
302 int rinfo = 0;
303 double rcond = ct_dgbcon('1', m_n, m_kl, m_ku, ludata.data(),
304 ldab, m_ipiv->data.data(), a1norm, work_.data(), iwork_.data(), rinfo);
305 if (rinfo != 0) {
306 throw CanteraError("BandMatrix::rcond", "DGBCON returned INFO = {}", rinfo);
307 }
308 return rcond;
309#else
310 throw CanteraError("BandMatrix::rcond", "not implemented when LAPACK is missing");
311#endif
312}
313
315{
316 return 0;
317}
318
320{
321 double value = 0.0;
322 for (size_t j = 0; j < m_n; j++) {
323 double sum = 0.0;
324 size_t start = (j >= m_ku) ? j - m_ku : 0;
325 size_t stop = std::min(j + m_kl + 1, m_n);
326 for (size_t i = start; i < stop; i++) {
327 sum += std::abs(_value(i,j));
328 }
329 value = std::max(sum, value);
330 }
331 return value;
332}
333
334size_t BandMatrix::checkRows(double& valueSmall) const
335{
336 valueSmall = 1.0E300;
337 size_t iSmall = npos;
338 for (size_t i = 0; i < m_n; i++) {
339 double valueS = 0.0;
340 size_t start = (i > m_kl) ? i - m_kl : 0;
341 size_t stop = std::min(i + m_ku + 1, m_n);
342 for (size_t j = start; j < stop; j++) {
343 valueS = std::max(fabs(_value(i,j)), valueS);
344 }
345 if (valueS < valueSmall) {
346 iSmall = i;
347 valueSmall = valueS;
348 if (valueSmall == 0.0) {
349 return iSmall;
350 }
351 }
352 }
353 return iSmall;
354}
355
356size_t BandMatrix::checkColumns(double& valueSmall) const
357{
358 valueSmall = 1.0E300;
359 size_t jSmall = npos;
360 for (size_t j = 0; j < m_n; j++) {
361 double valueS = 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 valueS = std::max(fabs(_value(i,j)), valueS);
366 }
367 if (valueS < valueSmall) {
368 jSmall = j;
369 valueSmall = valueS;
370 if (valueSmall == 0.0) {
371 return jSmall;
372 }
373 }
374 }
375 return jSmall;
376}
377
378}
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
double & operator()(size_t i, size_t j) override
Index into the (i,j) element.
vector< double > ludata
Factorized data.
Definition BandMatrix.h:204
size_t m_ku
Number of super diagonals of the matrix.
Definition BandMatrix.h:215
size_t m_kl
Number of subdiagonals of the matrix.
Definition BandMatrix.h:212
void mult(span< const double > b, span< double > prod) const override
Multiply A*b and write result to prod.
double m_zero
value of zero
Definition BandMatrix.h:218
double rcond(double a1norm) override
Returns an estimate of the inverse of the condition number for the matrix.
size_t nSubDiagonals() const
Number of subdiagonals.
void resize(size_t n, size_t kl, size_t ku, double v=0.0)
Resize the matrix problem.
size_t ldim() const
Return the number of rows of storage needed for the band storage.
void leftMult(span< const double > b, span< double > prod) const override
Multiply b*A and write result to prod.
void solve(span< const double > b, span< double > x)
Solve the matrix problem Ax = b.
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.
void bfill(double v=0.0)
Fill or zero the matrix.
BandMatrix()
Base Constructor.
size_t checkColumns(double &valueSmall) const override
Check to see if we have any zero columns in the Jacobian.
vector< double > work_
Extra dp work array needed - size = 3n.
Definition BandMatrix.h:229
unique_ptr< PivData > m_ipiv
Pivot vector.
Definition BandMatrix.h:223
size_t nColumns() const
Number of columns.
size_t checkRows(double &valueSmall) const override
Check to see if we have any zero rows in the Jacobian.
size_t nSuperDiagonals() const
Number of superdiagonals.
vector< double * > m_lu_col_ptrs
Pointers into ludata used by the SUNDIALS solver.
Definition BandMatrix.h:206
size_t nRows() const override
Return the number of rows in the matrix.
vector< int > iwork_
Extra work array needed - size = n.
Definition BandMatrix.h:226
vector< double > data
Matrix data.
Definition BandMatrix.h:201
void factor() override
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
void zero() override
Zero the matrix elements.
size_t m_n
Number of rows and columns of the matrix.
Definition BandMatrix.h:209
double oneNorm() const override
Returns the one norm of the matrix.
double _value(size_t i, size_t j) const
Return the value of the (i,j) element for (i,j) within the bandwidth.
double & value(size_t i, size_t j)
Return a changeable reference to element (i,j).
int factorAlgorithm() const override
Returns the factor algorithm used.
Base class for exceptions thrown by Cantera classes.
int m_factored
Indicates whether the matrix is factored.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:183
std::ostream & operator<<(std::ostream &s, const Array2D &m)
Output the current contents of the Array2D object.
Definition Array.cpp:101
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...