Cantera  3.3.0a1
Loading...
Searching...
No Matches
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
20DenseMatrix::DenseMatrix(size_t n, size_t m, double v) :
21 Array2D(n, m, v)
22{
23 m_ipiv.resize(std::max(n, m));
24 m_colPts.resize(m);
25 if (!m_data.empty()) {
26 for (size_t j = 0; j < m; j++) {
27 m_colPts[j] = &m_data[m_nrows*j];
28 }
29 }
30}
31
33 Array2D(y)
34{
35 m_ipiv = y.ipiv();
36 m_colPts.resize(m_ncols);
37 if (!m_data.empty()) {
38 for (size_t j = 0; j < m_ncols; j++) {
39 m_colPts[j] = &m_data[m_nrows*j];
40 }
41 }
42}
43
44DenseMatrix& DenseMatrix::operator=(const DenseMatrix& y)
45{
46 if (&y == this) {
47 return *this;
48 }
49 Array2D::operator=(y);
50 m_ipiv = y.ipiv();
51 m_colPts.resize(m_ncols);
52 for (size_t j = 0; j < m_ncols; j++) {
53 m_colPts[j] = &m_data[m_nrows*j];
54 }
55 return *this;
56}
57
58void DenseMatrix::resize(size_t n, size_t m, double v)
59{
60 Array2D::resize(n,m,v);
61 m_ipiv.resize(std::max(n,m));
62 m_colPts.resize(m_ncols);
63 if (!m_data.empty()) {
64 for (size_t j = 0; j < m_ncols; j++) {
65 m_colPts[j] = &m_data[m_nrows*j];
66 }
67 }
68}
69
70double* const* DenseMatrix::colPts()
71{
72 return &m_colPts[0];
73}
74
75const double* const* DenseMatrix::const_colPts() const
76{
77 return &m_colPts[0];
78}
79
80void DenseMatrix::mult(const double* b, double* prod) const
81{
82#if CT_USE_LAPACK
83 ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
84 static_cast<int>(nRows()),
85 static_cast<int>(nColumns()), 1.0, ptrColumn(0),
86 static_cast<int>(nRows()), b, 1, 0.0, prod, 1);
87#else
88 MappedMatrix mat(const_cast<double*>(m_data.data()), nRows(), nColumns());
89 MappedVector bm(const_cast<double*>(b), nColumns());
90 MappedVector pm(prod, nRows());
91 pm = mat * bm;
92#endif
93}
94
95void DenseMatrix::mult(const DenseMatrix& B, DenseMatrix& prod) const
96{
97 if (nColumns() != B.nRows()) {
98 throw CanteraError("DenseMatrix::mult",
99 "Inner matrix dimensions do not agree: {} != {}", nColumns(), B.nRows());
100 }
101 if (nRows() != prod.nRows() || B.nColumns() != prod.nColumns()) {
102 throw CanteraError("DenseMatrix::mult",
103 "Output matrix has wrong dimensions: {}x{} != {}x{}",
104 prod.nRows(), prod.nColumns(), nRows(), B.nColumns());
105 }
106 const double* const* bcols = B.const_colPts();
107 double* const* prodcols = prod.colPts();
108 for (size_t col=0; col < B.nColumns(); ++col) {
109 // Loop over ncols multiplying A*column of B and storing in
110 // corresponding prod column
111 mult(bcols[col], prodcols[col]);
112 }
113}
114
115void DenseMatrix::leftMult(const double* const b, double* const prod) const
116{
117 for (size_t n = 0; n < nColumns(); n++) {
118 double sum = 0.0;
119 for (size_t i = 0; i < nRows(); i++) {
120 sum += value(i,n)*b[i];
121 }
122 prod[n] = sum;
123 }
124}
125
126vector<int>& DenseMatrix::ipiv()
127{
128 return m_ipiv;
129}
130
131void solve(DenseMatrix& A, double* b, size_t nrhs, size_t ldb)
132{
133 if (A.nColumns() != A.nRows()) {
134 throw CanteraError("solve(DenseMatrix& A, double* b)", "Can only solve a square matrix");
135 }
136
137 int info = 0;
138 if (ldb == 0) {
139 ldb = A.nColumns();
140 }
141 #if CT_USE_LAPACK
142 ct_dgetrf(A.nRows(), A.nColumns(), A.ptrColumn(0),
143 A.nRows(), &A.ipiv()[0], info);
144 if (info > 0) {
145 throw CanteraError("solve(DenseMatrix& A, double* b)",
146 "DGETRF returned INFO = {}. U(i,i) is exactly zero. The factorization has"
147 " been completed, but the factor U is exactly singular, and division by zero will occur if "
148 "it is used to solve a system of equations.", info);
149 } else if (info < 0) {
150 throw CanteraError("solve(DenseMatrix& A, double* b)",
151 "DGETRF returned INFO = {}. The argument i has an illegal value", info);
152 }
153
154 ct_dgetrs(ctlapack::NoTranspose, A.nRows(), nrhs, A.ptrColumn(0),
155 A.nRows(), &A.ipiv()[0], b, ldb, info);
156 if (info != 0) {
157 throw CanteraError("solve(DenseMatrix& A, double* b)",
158 "DGETRS returned INFO = {}", info);
159 }
160 #else
161 MappedMatrix Am(&A(0,0), A.nRows(), A.nColumns());
162 #ifdef NDEBUG
163 auto lu = Am.partialPivLu();
164 #else
165 auto lu = Am.fullPivLu();
166 if (lu.nonzeroPivots() < static_cast<long int>(A.nColumns())) {
167 throw CanteraError("solve(DenseMatrix& A, double* b)",
168 "Matrix appears to be rank-deficient: non-zero pivots = {}; columns = {}",
169 lu.nonzeroPivots(), A.nColumns());
170 }
171 #endif
172 for (size_t i = 0; i < nrhs; i++) {
173 MappedVector bm(b + ldb*i, A.nColumns());
174 bm = lu.solve(bm);
175 }
176 #endif
177}
178
180{
181 solve(A, b.ptrColumn(0), b.nColumns(), b.nRows());
182}
183
184void multiply(const DenseMatrix& A, const double* const b, double* const prod)
185{
186 A.mult(b, prod);
187}
188
189void increment(const DenseMatrix& A, const double* b, double* prod)
190{
191 #if CT_USE_LAPACK
192 ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
193 static_cast<int>(A.nRows()), static_cast<int>(A.nColumns()), 1.0,
194 A.ptrColumn(0), static_cast<int>(A.nRows()), b, 1, 1.0, prod, 1);
195 #else
196 MappedMatrix Am(&const_cast<DenseMatrix&>(A)(0,0), A.nRows(), A.nColumns());
197 MappedVector bm(const_cast<double*>(b), A.nColumns());
198 MappedVector pm(prod, A.nRows());
199 pm += Am * bm;
200 #endif
201}
202
203void invert(DenseMatrix& A, size_t nn)
204{
205 int info=0;
206 #if CT_USE_LAPACK
207 integer n = static_cast<int>(nn != npos ? nn : A.nRows());
208 ct_dgetrf(n, n, A.ptrColumn(0), static_cast<int>(A.nRows()),
209 &A.ipiv()[0], info);
210 if (info != 0) {
211 throw CanteraError("invert(DenseMatrix& A, size_t nn)", "DGETRS returned INFO = {}", info);
212 }
213
214 vector<double> work(n);
215 integer lwork = static_cast<int>(work.size());
216 ct_dgetri(n, A.ptrColumn(0), static_cast<int>(A.nRows()),
217 &A.ipiv()[0], &work[0], lwork, info);
218 if (info != 0) {
219 throw CanteraError("invert(DenseMatrix& A, size_t nn)", "DGETRI returned INFO={}", info);
220 }
221 #else
222 MappedMatrix Am(&A(0,0), A.nRows(), A.nColumns());
223 if (nn == npos) {
224 Am = Am.inverse();
225 } else {
226 Am.topLeftCorner(nn, nn) = Am.topLeftCorner(nn, nn).inverse();
227 }
228 #endif
229}
230
231}
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:32
vector< double > m_data
Data stored in a single array.
Definition Array.h:219
size_t m_nrows
Number of rows.
Definition Array.h:222
size_t nRows() const
Number of rows.
Definition Array.h:176
size_t m_ncols
Number of columns.
Definition Array.h:225
size_t nColumns() const
Number of columns.
Definition Array.h:181
double * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are contiguous in memory.
Definition Array.h:203
double & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
Definition Array.h:160
virtual 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:47
Base class for exceptions thrown by Cantera classes.
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition DenseMatrix.h:42
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, double v=0.0) override
Resize the matrix.
vector< int > & ipiv()
Return a changeable value of the pivot vector.
DenseMatrix()=default
Default Constructor.
const double *const * const_colPts() const
Return a const vector of const pointers to the columns.
vector< double * > m_colPts
Vector of column pointers.
vector< int > m_ipiv
Vector of pivots. Length is equal to the max of m and n.
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
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.
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
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.
void solve(DenseMatrix &A, double *b, size_t nrhs, size_t ldb)
Solve Ax = b. Array b is overwritten on exit with x.
void invert(DenseMatrix &A, size_t nn)
invert A. A is overwritten with A^-1.
Contains declarations for string manipulation functions within Cantera.