13 #if CT_SUNDIALS_VERSION >= 30
14 #include "sunlinsol/sunlinsol_band.h"
16 #include "cvodes/cvodes_dense.h"
17 #include "cvodes/cvodes_band.h"
31struct BandMatrix::PivData {
34#elif CT_SUNDIALS_VERSION >= 30
35 std::vector<sunindextype>
data;
36#elif CT_SUNDIALS_VERSION >= 25
37 std::vector<long int>
data;
39 std::vector<int>
data;
48 m_ipiv{new PivData()},
53BandMatrix::~BandMatrix()
63 m_ipiv{new PivData()},
66 data.resize(n*(2*kl + ku + 1));
67 ludata.resize(n*(2*kl + ku + 1));
72 m_lu_col_ptrs.resize(n);
73 size_t ldab = (2*kl + ku + 1);
74 for (
size_t j = 0; j < n; j++) {
76 m_lu_col_ptrs[j] = &
ludata[ldab * j];
86 m_ipiv{new PivData()},
96 m_lu_col_ptrs.resize(
m_n);
98 for (
size_t j = 0; j <
m_n; j++) {
100 m_lu_col_ptrs[j] = &
ludata[ldab * j];
104BandMatrix& BandMatrix::operator=(
const BandMatrix& y)
109 GeneralMatrix::operator=(y);
113 m_ipiv->data = y.m_ipiv->data;
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++) {
121 m_lu_col_ptrs[j] = &
ludata[ldab * j];
132 data.resize(n*(2*kl + ku + 1));
133 ludata.resize(n*(2*kl + ku + 1));
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++) {
141 m_lu_col_ptrs[j] = &
ludata[ldab * j];
148 std::fill(
data.begin(),
data.end(), v);
154 std::fill(
data.begin(),
data.end(), 0.0);
171 if (i + m_ku < j || i > j +
m_kl) {
179 if (i + m_ku < j || i > j +
m_kl) {
222 for (
size_t m = 0; m <
m_n; m++) {
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];
235 for (
size_t n = 0; n <
m_n; n++) {
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];
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());
261 "Factorization failed with DGBTRF error code {}.", m_info);
284 m_ipiv->data.data(), b, ldb, m_info);
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);
296 "Linear solve failed with DGBTRS error code {}.", m_info);
325 for (
size_t i = 0; i < m.
nRows(); i++) {
327 for (
size_t j = 1; j < m.nColumns(); j++) {
341 throw CanteraError(
"BandMatrix::rcond",
"matrix isn't factored correctly");
350 throw CanteraError(
"BandMatrix::rcond",
"DGBCON returned INFO = {}", rinfo);
354 throw CanteraError(
"BandMatrix::rcond",
"not implemented when LAPACK is missing");
366 for (
size_t j = 0; j <
m_n; j++) {
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));
380 valueSmall = 1.0E300;
381 size_t iSmall =
npos;
382 for (
size_t i = 0; i <
m_n; i++) {
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);
389 if (valueS < valueSmall) {
392 if (valueSmall == 0.0) {
402 valueSmall = 1.0E300;
403 size_t jSmall =
npos;
404 for (
size_t j = 0; j <
m_n; j++) {
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);
411 if (valueS < valueSmall) {
414 if (valueSmall == 0.0) {
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.
virtual void zero()
Zero the matrix elements.
virtual doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j.
doublereal _value(size_t i, size_t j) const
Return the value of the (i,j) element for (i,j) within the bandwidth.
size_t m_ku
Number of super diagonals of the matrix.
size_t m_kl
Number of subdiagonals of the matrix.
vector_fp ludata
Factorized data.
virtual size_t checkRows(doublereal &valueSmall) const
Check to see if we have any zero rows in the Jacobian.
std::vector< doublereal * > m_colPtrs
Vector of column pointers.
size_t nSubDiagonals() const
Number of subdiagonals.
size_t ldim() const
Return the number of rows of storage needed for the band storage.
int factor()
Perform an LU decomposition, the LAPACK routine DGBTRF is used.
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.
virtual void mult(const doublereal *b, doublereal *prod) const
Multiply A*b and write result to prod.
virtual size_t nRows() const
Return the number of rows in the matrix.
virtual doublereal rcond(doublereal a1norm)
Returns an estimate of the inverse of the condition number for the matrix.
vector_int iwork_
Extra work array needed - size = n.
virtual vector_fp::iterator begin()
Returns an iterator for the start of the band storage data.
std::unique_ptr< PivData > m_ipiv
Pivot vector.
void bfill(doublereal v=0.0)
Fill or zero the matrix.
BandMatrix()
Base Constructor.
size_t nColumns() const
Number of columns.
doublereal m_zero
value of zero
doublereal & operator()(size_t i, size_t j)
Index into the (i,j) element.
size_t nSuperDiagonals() const
Number of superdiagonals.
vector_fp work_
Extra dp work array needed - size = 3n.
int solve(const doublereal *const b, doublereal *const x)
Solve the matrix problem Ax = b.
virtual int factorAlgorithm() const
Returns the factor algorithm used.
vector_fp::iterator end()
Returns an iterator for the end of the band storage data.
void resize(size_t n, size_t kl, size_t ku, doublereal v=0.0)
Resize the matrix problem.
virtual size_t checkColumns(doublereal &valueSmall) const
Check to see if we have any zero columns in the Jacobian.
virtual void leftMult(const doublereal *const b, doublereal *const prod) const
Multiply b*A and write result to prod.
virtual doublereal *const * colPts()
Return a vector of const pointers to the columns.
doublereal & value(size_t i, size_t j)
Return a changeable reference to element (i,j).
size_t m_n
Number of rows and columns of the matrix.
vector_fp data
Matrix data.
virtual doublereal oneNorm() const
Returns the one norm of the matrix.
Base class for exceptions thrown by Cantera classes.
int m_factored
Indicates whether the matrix is factored.
Namespace for the Cantera kernel.
const size_t npos
index returned by functions to indicate "no position"
std::vector< int > vector_int
Vector of ints.
std::ostream & operator<<(std::ostream &s, const Array2D &m)
Output the current contents of the Array2D object.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector operations (see Templated Utility Functions)...