13 #include "cantera/numerics/eigen_dense.h"
20 m_useReturnErrorCode(0),
27 m_useReturnErrorCode(0),
30 m_ipiv.resize(std::max(n, m));
33 for (
size_t j = 0; j < m; j++) {
41 m_useReturnErrorCode(0),
47 for (
size_t j = 0; j <
m_ncols; j++) {
53 DenseMatrix& DenseMatrix::operator=(
const DenseMatrix& y)
58 Array2D::operator=(y);
61 for (
size_t j = 0; j <
m_ncols; j++) {
72 m_ipiv.resize(std::max(n,m));
75 for (
size_t j = 0; j <
m_ncols; j++) {
81 doublereal*
const* DenseMatrix::colPts()
91 void DenseMatrix::mult(
const double* b,
double* prod)
const
94 ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
95 static_cast<int>(
nRows()),
97 static_cast<int>(
nRows()), b, 1, 0.0, prod, 1);
100 MappedVector bm(
const_cast<double*
>(b),
nColumns());
101 MappedVector pm(prod,
nRows());
110 "Inner matrix dimensions do not agree: {} != {}",
nColumns(), B.
nRows());
114 "Output matrix has wrong dimensions: {}x{} != {}x{}",
118 doublereal*
const* prodcols = prod.colPts();
119 for (
size_t col=0; col < B.
nColumns(); ++col) {
122 mult(bcols[col], prodcols[col]);
128 for (
size_t n = 0; n <
nColumns(); n++) {
130 for (
size_t i = 0; i <
nRows(); i++) {
131 sum +=
value(i,n)*b[i];
146 writelogf(
"solve(DenseMatrix& A, double* b): Can only solve a square matrix\n");
148 throw CanteraError(
"solve(DenseMatrix& A, double* b)",
"Can only solve a square matrix");
160 writelogf(
"solve(DenseMatrix& A, double* b): DGETRF returned INFO = %d U(i,i) is exactly zero. The factorization has"
161 " been completed, but the factor U is exactly singular, and division by zero will occur if "
162 "it is used to solve a system of equations.\n", info);
166 "DGETRF returned INFO = {}. U(i,i) is exactly zero. The factorization has"
167 " been completed, but the factor U is exactly singular, and division by zero will occur if "
168 "it is used to solve a system of equations.", info);
171 }
else if (info < 0) {
173 writelogf(
"solve(DenseMatrix& A, double* b): DGETRF returned INFO = %d. The argument i has an illegal value\n", info);
177 "DGETRF returned INFO = {}. The argument i has an illegal value", info);
184 writelog(
"solve(DenseMatrix& A, double* b): DGETRS returned INFO = {}\n", info);
188 "DGETRS returned INFO = {}", info);
194 auto lu = Am.partialPivLu();
196 auto lu = Am.fullPivLu();
197 if (lu.nonzeroPivots() <
static_cast<long int>(A.
nColumns())) {
199 "Matrix appears to be rank-deficient: non-zero pivots = {}; columns = {}",
203 for (
size_t i = 0; i < nrhs; i++) {
204 MappedVector bm(b + ldb*i, A.
nColumns());
224 ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
225 static_cast<int>(A.
nRows()),
static_cast<int>(A.
nColumns()), 1.0,
229 MappedVector bm(
const_cast<double*
>(b), A.
nColumns());
230 MappedVector pm(prod, A.
nRows());
239 integer n =
static_cast<int>(nn !=
npos ? nn : A.
nRows());
244 writelogf(
"invert(DenseMatrix& A, size_t nn): DGETRS returned INFO = %d\n", info);
247 throw CanteraError(
"invert(DenseMatrix& A, size_t nn)",
"DGETRS returned INFO = {}", info);
253 integer lwork =
static_cast<int>(work.size());
255 &A.
ipiv()[0], &work[0], lwork, info);
258 writelogf(
"invert(DenseMatrix& A, size_t nn): DGETRS returned INFO = %d\n", info);
261 throw CanteraError(
"invert(DenseMatrix& A, size_t nn)",
"DGETRI returned INFO={}", info);
269 Am.topLeftCorner(nn, nn) = Am.topLeftCorner(nn, nn).inverse();
Headers for the DenseMatrix object, which deals with dense rectangular matrices and description of th...
A class for 2D arrays stored in column-major (Fortran-compatible) form.
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the array, and fill the new entries with 'v'.
size_t m_nrows
Number of rows.
vector_fp m_data
Data stored in a single array.
doublereal & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
size_t nRows() const
Number of rows.
size_t m_ncols
Number of columns.
doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are contiguous in memory.
size_t nColumns() const
Number of columns.
Base class for exceptions thrown by Cantera classes.
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
virtual void leftMult(const double *const b, double *const prod) const
Left-multiply the matrix by transpose(b), and write the result to prod.
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the matrix.
std::vector< doublereal * > m_colPts
Vector of column pointers.
int m_printLevel
Print Level.
int m_useReturnErrorCode
Error Handling Flag.
vector_int m_ipiv
Vector of pivots. Length is equal to the max of m and n.
const doublereal *const * const_colPts() const
Return a const vector of const pointers to the columns.
DenseMatrix()
Default Constructor.
vector_int & ipiv()
Return a changeable value of the pivot vector.
const size_t npos
index returned by functions to indicate "no position"
std::vector< int > vector_int
Vector of ints.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
void writelog(const std::string &fmt, const Args &... args)
Write a formatted message to the screen.
void writelogf(const char *fmt, const Args &... args)
Write a formatted message to the screen.
Namespace for the Cantera kernel.
void increment(const DenseMatrix &A, const double *b, double *prod)
Multiply A*b and add it to the result in prod. Uses BLAS routine DGEMV.
void multiply(const DenseMatrix &A, const double *const b, double *const prod)
Multiply A*b and return the result in prod. Uses BLAS routine DGEMV.
int invert(DenseMatrix &A, size_t nn)
invert A. A is overwritten with A^-1.
int solve(DenseMatrix &A, double *b, size_t nrhs, size_t ldb)
Solve Ax = b. Array b is overwritten on exit with x.
Contains declarations for string manipulation functions within Cantera.