13 #include "sunlinsol/sunlinsol_band.h"
30 vector<sunindextype> data;
39BandMatrix::~BandMatrix()
50 data.resize(n*(2*kl + ku + 1));
51 ludata.resize(n*(2*kl + ku + 1));
56 size_t ldab = (2*kl + ku + 1);
57 for (
size_t j = 0; j < n; j++) {
69 m_ipiv{new PivData()},
75 for (
size_t j = 0; j <
m_n; j++) {
80BandMatrix& BandMatrix::operator=(
const BandMatrix& y)
85 GeneralMatrix::operator=(y);
89 m_ipiv->data = y.m_ipiv->data;
94 for (
size_t j = 0; j <
m_n; j++) {
106 data.resize(n*(2*kl + ku + 1));
107 ludata.resize(n*(2*kl + ku + 1));
110 for (
size_t j = 0; j <
m_n; j++) {
120 std::fill(
data.begin(),
data.end(), v);
126 std::fill(
data.begin(),
data.end(), 0.0);
143 if (i + m_ku < j || i > j +
m_kl) {
151 if (i + m_ku < j || i > j +
m_kl) {
196 for (
size_t m = 0; m <
m_n; m++) {
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];
211 for (
size_t n = 0; n <
m_n; n++) {
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];
231 long int smu = nu + nl;
234 nu, nl, smu,
m_ipiv->data.data());
238 "Factorization failed with DGBTRF error code {}.", m_info);
247 copy(b.begin(), b.begin() +
m_n, x.begin());
263 m_ipiv->data.data(), b.data(), ldb, m_info);
266 "Linear solve failed with DGBTRS error code {}.", m_info);
271 long int smu = nu + nl;
273 SUNDlsMat_bandGBTRS(a,
static_cast<long int>(
nColumns()), smu, nl,
274 m_ipiv->data.data(), b.data());
281 for (
size_t i = 0; i < m.
nRows(); i++) {
283 for (
size_t j = 1; j < m.nColumns(); j++) {
297 throw CanteraError(
"BandMatrix::rcond",
"matrix isn't factored correctly");
306 throw CanteraError(
"BandMatrix::rcond",
"DGBCON returned INFO = {}", rinfo);
310 throw CanteraError(
"BandMatrix::rcond",
"not implemented when LAPACK is missing");
322 for (
size_t j = 0; j <
m_n; j++) {
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));
336 valueSmall = 1.0E300;
337 size_t iSmall =
npos;
338 for (
size_t i = 0; i <
m_n; i++) {
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);
345 if (valueS < valueSmall) {
348 if (valueSmall == 0.0) {
358 valueSmall = 1.0E300;
359 size_t jSmall =
npos;
360 for (
size_t j = 0; j <
m_n; j++) {
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);
367 if (valueS < valueSmall) {
370 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.
double & operator()(size_t i, size_t j) override
Index into the (i,j) element.
vector< double > ludata
Factorized data.
size_t m_ku
Number of super diagonals of the matrix.
size_t m_kl
Number of subdiagonals of the matrix.
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
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.
unique_ptr< PivData > m_ipiv
Pivot vector.
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.
size_t nRows() const override
Return the number of rows in the matrix.
vector< int > iwork_
Extra work array needed - size = n.
vector< double > data
Matrix data.
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.
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.
const size_t npos
index returned by functions to indicate "no position"
std::ostream & operator<<(std::ostream &s, const Array2D &m)
Output the current contents of the Array2D object.
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...