Cantera 2.6.0
DenseMatrix.cpp
Go to the documentation of this file.
1/**
2 * @file DenseMatrix.cpp
3 */
4
5// This file is part of Cantera. See License.txt in the top-level directory or
6// at https://cantera.org/license.txt for license and copyright information.
7
10#include "cantera/base/global.h"
11#if CT_USE_LAPACK
13#else
14 #include "cantera/numerics/eigen_dense.h"
15#endif
16
17namespace Cantera
18{
19
21 m_useReturnErrorCode(0),
22 m_printLevel(0)
23{
24}
25
26DenseMatrix::DenseMatrix(size_t n, size_t m, doublereal v) :
27 Array2D(n, m, v),
28 m_useReturnErrorCode(0),
29 m_printLevel(0)
30{
31 m_ipiv.resize(std::max(n, m));
32 m_colPts.resize(m);
33 if (!m_data.empty()) {
34 for (size_t j = 0; j < m; j++) {
35 m_colPts[j] = &m_data[m_nrows*j];
36 }
37 }
38}
39
41 Array2D(y),
42 m_useReturnErrorCode(0),
43 m_printLevel(0)
44{
45 m_ipiv = y.ipiv();
46 m_colPts.resize(m_ncols);
47 if (!m_data.empty()) {
48 for (size_t j = 0; j < m_ncols; j++) {
49 m_colPts[j] = &m_data[m_nrows*j];
50 }
51 }
52}
53
54DenseMatrix& DenseMatrix::operator=(const DenseMatrix& y)
55{
56 if (&y == this) {
57 return *this;
58 }
59 Array2D::operator=(y);
60 m_ipiv = y.ipiv();
61 m_colPts.resize(m_ncols);
62 for (size_t j = 0; j < m_ncols; j++) {
63 m_colPts[j] = &m_data[m_nrows*j];
64 }
65 m_useReturnErrorCode = y.m_useReturnErrorCode;
66 m_printLevel = y.m_printLevel;
67 return *this;
68}
69
70void DenseMatrix::resize(size_t n, size_t m, doublereal v)
71{
72 Array2D::resize(n,m,v);
73 m_ipiv.resize(std::max(n,m));
74 m_colPts.resize(m_ncols);
75 if (!m_data.empty()) {
76 for (size_t j = 0; j < m_ncols; j++) {
77 m_colPts[j] = &m_data[m_nrows*j];
78 }
79 }
80}
81
82doublereal* const* DenseMatrix::colPts()
83{
84 return &m_colPts[0];
85}
86
87const doublereal* const* DenseMatrix::const_colPts() const
88{
89 return &m_colPts[0];
90}
91
92void DenseMatrix::mult(const double* b, double* prod) const
93{
94#if CT_USE_LAPACK
95 ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
96 static_cast<int>(nRows()),
97 static_cast<int>(nColumns()), 1.0, ptrColumn(0),
98 static_cast<int>(nRows()), b, 1, 0.0, prod, 1);
99#else
100 MappedMatrix mat(const_cast<double*>(m_data.data()), nRows(), nColumns());
101 MappedVector bm(const_cast<double*>(b), nColumns());
102 MappedVector pm(prod, nRows());
103 pm = mat * bm;
104#endif
105}
106
107void DenseMatrix::mult(const DenseMatrix& B, DenseMatrix& prod) const
108{
109 if (nColumns() != B.nRows()) {
110 throw CanteraError("DenseMatrix::mult",
111 "Inner matrix dimensions do not agree: {} != {}", nColumns(), B.nRows());
112 }
113 if (nRows() != prod.nRows() || B.nColumns() != prod.nColumns()) {
114 throw CanteraError("DenseMatrix::mult",
115 "Output matrix has wrong dimensions: {}x{} != {}x{}",
116 prod.nRows(), prod.nColumns(), nRows(), B.nColumns());
117 }
118 const doublereal* const* bcols = B.const_colPts();
119 doublereal* const* prodcols = prod.colPts();
120 for (size_t col=0; col < B.nColumns(); ++col) {
121 // Loop over ncols multiplying A*column of B and storing in
122 // corresponding prod column
123 mult(bcols[col], prodcols[col]);
124 }
125}
126
127void DenseMatrix::leftMult(const double* const b, double* const prod) const
128{
129 for (size_t n = 0; n < nColumns(); n++) {
130 double sum = 0.0;
131 for (size_t i = 0; i < nRows(); i++) {
132 sum += value(i,n)*b[i];
133 }
134 prod[n] = sum;
135 }
136}
137
139{
140 return m_ipiv;
141}
142
143int solve(DenseMatrix& A, double* b, size_t nrhs, size_t ldb)
144{
145 if (A.nColumns() != A.nRows()) {
146 if (A.m_printLevel) {
147 writelogf("solve(DenseMatrix& A, double* b): Can only solve a square matrix\n");
148 }
149 throw CanteraError("solve(DenseMatrix& A, double* b)", "Can only solve a square matrix");
150 }
151
152 int info = 0;
153 if (ldb == 0) {
154 ldb = A.nColumns();
155 }
156 #if CT_USE_LAPACK
157 ct_dgetrf(A.nRows(), A.nColumns(), A.ptrColumn(0),
158 A.nRows(), &A.ipiv()[0], info);
159 if (info > 0) {
160 if (A.m_printLevel) {
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);
164 }
165 if (!A.m_useReturnErrorCode) {
166 throw CanteraError("solve(DenseMatrix& A, double* b)",
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);
170 }
171 return info;
172 } else if (info < 0) {
173 if (A.m_printLevel) {
174 writelogf("solve(DenseMatrix& A, double* b): DGETRF returned INFO = %d. The argument i has an illegal value\n", info);
175 }
176
177 throw CanteraError("solve(DenseMatrix& A, double* b)",
178 "DGETRF returned INFO = {}. The argument i has an illegal value", info);
179 }
180
181 ct_dgetrs(ctlapack::NoTranspose, A.nRows(), nrhs, A.ptrColumn(0),
182 A.nRows(), &A.ipiv()[0], b, ldb, info);
183 if (info != 0) {
184 if (A.m_printLevel) {
185 writelog("solve(DenseMatrix& A, double* b): DGETRS returned INFO = {}\n", info);
186 }
187 if (info < 0 || !A.m_useReturnErrorCode) {
188 throw CanteraError("solve(DenseMatrix& A, double* b)",
189 "DGETRS returned INFO = {}", info);
190 }
191 }
192 #else
193 MappedMatrix Am(&A(0,0), A.nRows(), A.nColumns());
194 #ifdef NDEBUG
195 auto lu = Am.partialPivLu();
196 #else
197 auto lu = Am.fullPivLu();
198 if (lu.nonzeroPivots() < static_cast<long int>(A.nColumns())) {
199 throw CanteraError("solve(DenseMatrix& A, double* b)",
200 "Matrix appears to be rank-deficient: non-zero pivots = {}; columns = {}",
201 lu.nonzeroPivots(), A.nColumns());
202 }
203 #endif
204 for (size_t i = 0; i < nrhs; i++) {
205 MappedVector bm(b + ldb*i, A.nColumns());
206 bm = lu.solve(bm);
207 }
208 #endif
209 return info;
210}
211
213{
214 return solve(A, b.ptrColumn(0), b.nColumns(), b.nRows());
215}
216
217void multiply(const DenseMatrix& A, const double* const b, double* const prod)
218{
219 A.mult(b, prod);
220}
221
222void increment(const DenseMatrix& A, const double* b, double* prod)
223{
224 #if CT_USE_LAPACK
225 ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
226 static_cast<int>(A.nRows()), static_cast<int>(A.nColumns()), 1.0,
227 A.ptrColumn(0), static_cast<int>(A.nRows()), b, 1, 1.0, prod, 1);
228 #else
229 MappedMatrix Am(&const_cast<DenseMatrix&>(A)(0,0), A.nRows(), A.nColumns());
230 MappedVector bm(const_cast<double*>(b), A.nColumns());
231 MappedVector pm(prod, A.nRows());
232 pm += Am * bm;
233 #endif
234}
235
236int invert(DenseMatrix& A, size_t nn)
237{
238 int info=0;
239 #if CT_USE_LAPACK
240 integer n = static_cast<int>(nn != npos ? nn : A.nRows());
241 ct_dgetrf(n, n, A.ptrColumn(0), static_cast<int>(A.nRows()),
242 &A.ipiv()[0], info);
243 if (info != 0) {
244 if (A.m_printLevel) {
245 writelogf("invert(DenseMatrix& A, size_t nn): DGETRS returned INFO = %d\n", info);
246 }
247 if (! A.m_useReturnErrorCode) {
248 throw CanteraError("invert(DenseMatrix& A, size_t nn)", "DGETRS returned INFO = {}", info);
249 }
250 return info;
251 }
252
253 vector_fp work(n);
254 integer lwork = static_cast<int>(work.size());
255 ct_dgetri(n, A.ptrColumn(0), static_cast<int>(A.nRows()),
256 &A.ipiv()[0], &work[0], lwork, info);
257 if (info != 0) {
258 if (A.m_printLevel) {
259 writelogf("invert(DenseMatrix& A, size_t nn): DGETRS returned INFO = %d\n", info);
260 }
261 if (! A.m_useReturnErrorCode) {
262 throw CanteraError("invert(DenseMatrix& A, size_t nn)", "DGETRI returned INFO={}", info);
263 }
264 }
265 #else
266 MappedMatrix Am(&A(0,0), A.nRows(), A.nColumns());
267 if (nn == npos) {
268 Am = Am.inverse();
269 } else {
270 Am.topLeftCorner(nn, nn) = Am.topLeftCorner(nn, nn).inverse();
271 }
272 #endif
273 return info;
274}
275
276}
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.
Definition: Array.h:30
size_t m_nrows
Number of rows.
Definition: Array.h:252
vector_fp m_data
Data stored in a single array.
Definition: Array.h:249
doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are contiguous in memory.
Definition: Array.h:233
size_t nRows() const
Number of rows.
Definition: Array.h:188
size_t m_ncols
Number of columns.
Definition: Array.h:255
size_t nColumns() const
Number of columns.
Definition: Array.h:193
doublereal & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
Definition: Array.h:172
void resize(size_t n, size_t m, double v=0.0)
Resize the array, and fill the new entries with 'v'.
Definition: Array.cpp:52
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition: DenseMatrix.h:51
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.
Definition: DenseMatrix.cpp:70
std::vector< doublereal * > m_colPts
Vector of column pointers.
Definition: DenseMatrix.h:128
int m_printLevel
Print Level.
Definition: DenseMatrix.h:149
int m_useReturnErrorCode
Error Handling Flag.
Definition: DenseMatrix.h:140
vector_int m_ipiv
Vector of pivots. Length is equal to the max of m and n.
Definition: DenseMatrix.h:125
const doublereal *const * const_colPts() const
Return a const vector of const pointers to the columns.
Definition: DenseMatrix.cpp:87
DenseMatrix()
Default Constructor.
Definition: DenseMatrix.cpp:20
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.
Definition: global.h:164
void writelogf(const char *fmt, const Args &... args)
Write a formatted message to the screen.
Definition: global.h:185
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
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"
Definition: ct_defs.h:192
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.
Definition: ct_defs.h:186
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:184
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.