14 #include "cantera/numerics/eigen_dense.h"
21 m_useReturnErrorCode(0),
28 m_useReturnErrorCode(0),
31 m_ipiv.resize(std::max(n, m));
34 for (
size_t j = 0; j < m; j++) {
42 m_useReturnErrorCode(0),
48 for (
size_t j = 0; j <
m_ncols; j++) {
54DenseMatrix& DenseMatrix::operator=(
const DenseMatrix& y)
59 Array2D::operator=(y);
62 for (
size_t j = 0; j <
m_ncols; j++) {
73 m_ipiv.resize(std::max(n,m));
76 for (
size_t j = 0; j <
m_ncols; j++) {
82doublereal*
const* DenseMatrix::colPts()
92void DenseMatrix::mult(
const double* b,
double* prod)
const
95 ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
96 static_cast<int>(
nRows()),
98 static_cast<int>(
nRows()), b, 1, 0.0, prod, 1);
101 MappedVector bm(
const_cast<double*
>(b),
nColumns());
102 MappedVector pm(prod,
nRows());
111 "Inner matrix dimensions do not agree: {} != {}",
nColumns(), B.
nRows());
115 "Output matrix has wrong dimensions: {}x{} != {}x{}",
119 doublereal*
const* prodcols = prod.colPts();
120 for (
size_t col=0; col < B.
nColumns(); ++col) {
123 mult(bcols[col], prodcols[col]);
129 for (
size_t n = 0; n <
nColumns(); n++) {
131 for (
size_t i = 0; i <
nRows(); i++) {
132 sum +=
value(i,n)*b[i];
147 writelogf(
"solve(DenseMatrix& A, double* b): Can only solve a square matrix\n");
149 throw CanteraError(
"solve(DenseMatrix& A, double* b)",
"Can only solve a square matrix");
161 writelogf(
"solve(DenseMatrix& A, double* b): DGETRF returned INFO = %d U(i,i) is exactly zero. The factorization has"
162 " been completed, but the factor U is exactly singular, and division by zero will occur if "
163 "it is used to solve a system of equations.\n", info);
167 "DGETRF returned INFO = {}. U(i,i) is exactly zero. The factorization has"
168 " been completed, but the factor U is exactly singular, and division by zero will occur if "
169 "it is used to solve a system of equations.", info);
172 }
else if (info < 0) {
174 writelogf(
"solve(DenseMatrix& A, double* b): DGETRF returned INFO = %d. The argument i has an illegal value\n", info);
178 "DGETRF returned INFO = {}. The argument i has an illegal value", info);
185 writelog(
"solve(DenseMatrix& A, double* b): DGETRS returned INFO = {}\n", info);
189 "DGETRS returned INFO = {}", info);
195 auto lu = Am.partialPivLu();
197 auto lu = Am.fullPivLu();
198 if (lu.nonzeroPivots() <
static_cast<long int>(A.
nColumns())) {
200 "Matrix appears to be rank-deficient: non-zero pivots = {}; columns = {}",
204 for (
size_t i = 0; i < nrhs; i++) {
205 MappedVector bm(b + ldb*i, A.
nColumns());
225 ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
226 static_cast<int>(A.
nRows()),
static_cast<int>(A.
nColumns()), 1.0,
230 MappedVector bm(
const_cast<double*
>(b), A.
nColumns());
231 MappedVector pm(prod, A.
nRows());
240 integer n =
static_cast<int>(nn !=
npos ? nn : A.
nRows());
245 writelogf(
"invert(DenseMatrix& A, size_t nn): DGETRS returned INFO = %d\n", info);
248 throw CanteraError(
"invert(DenseMatrix& A, size_t nn)",
"DGETRS returned INFO = {}", info);
254 integer lwork =
static_cast<int>(work.size());
256 &A.
ipiv()[0], &work[0], lwork, info);
259 writelogf(
"invert(DenseMatrix& A, size_t nn): DGETRS returned INFO = %d\n", info);
262 throw CanteraError(
"invert(DenseMatrix& A, size_t nn)",
"DGETRI returned INFO={}", info);
270 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.
size_t m_nrows
Number of rows.
vector_fp m_data
Data stored in a single array.
doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are contiguous in memory.
size_t nRows() const
Number of rows.
size_t m_ncols
Number of columns.
size_t nColumns() const
Number of columns.
doublereal & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
void resize(size_t n, size_t m, double v=0.0)
Resize the array, and fill the new entries with 'v'.
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.
This file contains definitions for utility functions and text for modules, inputfiles,...
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 *const b, double *const prod)
Multiply A*b and add it to the result in prod. Uses BLAS routine DGEMV.
const size_t npos
index returned by functions to indicate "no position"
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=npos)
invert A. A is overwritten with A^-1.
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.
int solve(DenseMatrix &A, double *b, size_t nrhs=1, size_t ldb=0)
Solve Ax = b. Array b is overwritten on exit with x.
Contains declarations for string manipulation functions within Cantera.