Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DenseMatrix.cpp
Go to the documentation of this file.
1 /**
2  * @file DenseMatrix.cpp
3  */
4 
5 // Copyright 2001 California Institute of Technology
6 
10 
11 namespace Cantera
12 {
13 
15  m_useReturnErrorCode(0),
16  m_printLevel(0)
17 {
18 }
19 
20 DenseMatrix::DenseMatrix(size_t n, size_t m, doublereal v) :
21  Array2D(n, m, v),
22  m_useReturnErrorCode(0),
23  m_printLevel(0)
24 {
25  m_ipiv.resize(std::max(n, m));
26  m_colPts.resize(m);
27  if (!m_data.empty()) {
28  for (size_t j = 0; j < m; j++) {
29  m_colPts[j] = &(m_data[m_nrows*j]);
30  }
31  }
32 }
33 
35  Array2D(y),
36  m_useReturnErrorCode(0),
37  m_printLevel(0)
38 {
39  m_ipiv = y.ipiv();
40  m_colPts.resize(m_ncols);
41  if (!m_data.empty()) {
42  for (size_t j = 0; j < m_ncols; j++) {
43  m_colPts[j] = &(m_data[m_nrows*j]);
44  }
45  }
46 }
47 
49 {
50  if (&y == this) {
51  return *this;
52  }
54  m_ipiv = y.ipiv();
55  m_colPts.resize(m_ncols);
56  for (size_t j = 0; j < m_ncols; j++) {
57  m_colPts[j] = &(m_data[m_nrows*j]);
58  }
61  return *this;
62 }
63 
64 void DenseMatrix::resize(size_t n, size_t m, doublereal v)
65 {
66  Array2D::resize(n,m,v);
67  m_ipiv.resize(std::max(n,m));
68  m_colPts.resize(m_ncols);
69  if (!m_data.empty()) {
70  for (size_t j = 0; j < m_ncols; j++) {
71  m_colPts[j] = &(m_data[m_nrows*j]);
72  }
73  }
74 }
75 
76 doublereal* const* DenseMatrix::colPts()
77 {
78  return &(m_colPts[0]);
79 }
80 
81 const doublereal* const* DenseMatrix::const_colPts() const
82 {
83  return &(m_colPts[0]);
84 }
85 
86 void DenseMatrix::mult(const double* b, double* prod) const
87 {
88  ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
89  static_cast<int>(nRows()),
90  static_cast<int>(nRows()), 1.0, ptrColumn(0),
91  static_cast<int>(nRows()), b, 1, 0.0, prod, 1);
92 }
93 
94 void DenseMatrix::mult(const DenseMatrix& B, DenseMatrix& prod) const
95 {
96  if (m_ncols != B.nColumns() || m_nrows != B.nRows() || m_ncols != m_nrows || m_ncols != prod.nColumns() || m_nrows != prod.nColumns()) {
97  throw CanteraError("mult(const DenseMatrix &B, DenseMatrix &prod)",
98  "Cannot multiply matrices that are not square and/or not the same size.");
99  }
100  const doublereal* const* bcols = B.const_colPts();
101  doublereal* const* prodcols = prod.colPts();
102  for (size_t col=0; col < m_ncols; ++col) {
103  // Loop over ncols multiplying A*column of B and storing in corresponding prod column
104  mult(bcols[col], prodcols[col]);
105  }
106 }
107 
108 void DenseMatrix::leftMult(const double* const b, double* const prod) const
109 {
110  size_t nc = nColumns();
111  size_t nr = nRows();
112  double sum = 0.0;
113  for (size_t n = 0; n < nc; n++) {
114  sum = 0.0;
115  for (size_t i = 0; i < nr; i++) {
116  sum += value(i,n)*b[i];
117  }
118  prod[n] = sum;
119  }
120 }
121 
123 {
124  return m_ipiv;
125 }
126 
127 int solve(DenseMatrix& A, double* b, size_t nrhs, size_t ldb)
128 {
129  int info = 0;
130  if (A.nColumns() != A.nRows()) {
131  if (A.m_printLevel) {
132  writelogf("solve(DenseMatrix& A, double* b): Can only solve a square matrix\n");
133  }
134  throw CELapackError("solve(DenseMatrix& A, double* b)", "Can only solve a square matrix");
135  }
136  ct_dgetrf(A.nRows(), A.nColumns(), A.ptrColumn(0),
137  A.nRows(), &A.ipiv()[0], info);
138  if (info > 0) {
139  if (A.m_printLevel) {
140  writelogf("solve(DenseMatrix& A, double* b): DGETRF returned INFO = %d U(i,i) is exactly zero. The factorization has"
141  " been completed, but the factor U is exactly singular, and division by zero will occur if "
142  "it is used to solve a system of equations.\n", info);
143  }
144  if (!A.m_useReturnErrorCode) {
145  throw CELapackError("solve(DenseMatrix& A, double* b)",
146  "DGETRF returned INFO = "+int2str(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.");
149  }
150  return info;
151  } else if (info < 0) {
152  if (A.m_printLevel) {
153  writelogf("solve(DenseMatrix& A, double* b): DGETRF returned INFO = %d. The argument i has an illegal value\n", info);
154  }
155 
156  throw CELapackError("solve(DenseMatrix& A, double* b)",
157  "DGETRF returned INFO = "+int2str(info) + ". The argument i has an illegal value");
158  }
159 
160  if (ldb == 0) {
161  ldb = A.nColumns();
162  }
163  ct_dgetrs(ctlapack::NoTranspose, A.nRows(), nrhs, A.ptrColumn(0),
164  A.nRows(), &A.ipiv()[0], b, ldb, info);
165  if (info != 0) {
166  if (A.m_printLevel) {
167  writelogf("solve(DenseMatrix& A, double* b): DGETRS returned INFO = %d\n", info);
168  }
169  if (info < 0 || !A.m_useReturnErrorCode) {
170  throw CELapackError("solve(DenseMatrix& A, double* b)", "DGETRS returned INFO = "+int2str(info));
171  }
172  }
173  return info;
174 }
175 
177 {
178  return solve(A, b.ptrColumn(0), b.nColumns(), b.nRows());
179 }
180 
181 void multiply(const DenseMatrix& A, const double* const b, double* const prod)
182 {
183  ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
184  static_cast<int>(A.nRows()), static_cast<int>(A.nColumns()), 1.0,
185  A.ptrColumn(0), static_cast<int>(A.nRows()), b, 1, 0.0, prod, 1);
186 }
187 
188 void increment(const DenseMatrix& A, const double* b, double* prod)
189 {
190  ct_dgemv(ctlapack::ColMajor, ctlapack::NoTranspose,
191  static_cast<int>(A.nRows()), static_cast<int>(A.nRows()), 1.0,
192  A.ptrColumn(0), static_cast<int>(A.nRows()), b, 1, 1.0, prod, 1);
193 }
194 
195 int invert(DenseMatrix& A, size_t nn)
196 {
197  integer n = static_cast<int>(nn != npos ? nn : A.nRows());
198  int info=0;
199  ct_dgetrf(n, n, A.ptrColumn(0), static_cast<int>(A.nRows()),
200  &A.ipiv()[0], info);
201  if (info != 0) {
202  if (A.m_printLevel) {
203  writelogf("invert(DenseMatrix& A, int nn): DGETRS returned INFO = %d\n", info);
204  }
205  if (! A.m_useReturnErrorCode) {
206  throw CELapackError("invert(DenseMatrix& A, int nn)", "DGETRS returned INFO = "+int2str(info));
207  }
208  return info;
209  }
210 
211  vector_fp work(n);
212  integer lwork = static_cast<int>(work.size());
213  ct_dgetri(n, A.ptrColumn(0), static_cast<int>(A.nRows()),
214  &A.ipiv()[0], &work[0], lwork, info);
215  if (info != 0) {
216  if (A.m_printLevel) {
217  writelogf("invert(DenseMatrix& A, int nn): DGETRS returned INFO = %d\n", info);
218  }
219  if (! A.m_useReturnErrorCode) {
220  throw CELapackError("invert(DenseMatrix& A, int nn)", "DGETRI returned INFO="+int2str(info));
221  }
222  }
223  return info;
224 }
225 
226 }
const doublereal *const * const_colPts() const
Return a const vector of const pointers to the columns.
Definition: DenseMatrix.cpp:81
DenseMatrix & operator=(const DenseMatrix &y)
Assignment operator.
Definition: DenseMatrix.cpp:48
size_t nRows() const
Number of rows.
Definition: Array.h:312
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:39
int invert(DenseMatrix &A, size_t nn)
invert A. A is overwritten with A^-1.
vector_fp m_data
Data stored in a single array.
Definition: Array.h:375
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the array, and fill the new entries with 'v'.
Definition: Array.h:121
Array2D & operator=(const Array2D &y)
assignment operator
Definition: Array.h:104
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:165
vector_int m_ipiv
Vector of pivots. Length is equal to the max of m and n.
Definition: DenseMatrix.h:155
int solve(DenseMatrix &A, double *b, size_t nrhs, size_t ldb)
Solve Ax = b. Array b is overwritten on exit with x.
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.
int m_printLevel
Print Level.
Definition: DenseMatrix.h:176
doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are contiguous in memory.
Definition: Array.h:358
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:29
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:159
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the matrix.
Definition: DenseMatrix.cpp:64
void writelogf(const char *fmt,...)
Write a formatted message to the screen.
Definition: global.cpp:38
std::vector< doublereal * > m_colPts
Vector of column pointers.
Definition: DenseMatrix.h:158
size_t nColumns() const
Number of columns.
Definition: Array.h:317
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.
size_t m_ncols
Number of columns.
Definition: Array.h:381
size_t m_nrows
Number of rows.
Definition: Array.h:378
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:99
doublereal & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
Definition: Array.h:295
virtual void leftMult(const double *const b, double *const prod) const
Left-multiply the matrix by transpose(b), and write the result to prod.
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:157
Headers for the DenseMatrix object, which deals with dense rectangular matrices and description of th...
Contains declarations for string manipulation functions within Cantera.
Exception thrown when an LAPACK error is encountered associated with inverting or solving a matrix...
Definition: DenseMatrix.h:35
DenseMatrix()
Default Constructor.
Definition: DenseMatrix.cpp:14
int m_useReturnErrorCode
Error Handling Flag.
Definition: DenseMatrix.h:168
vector_int & ipiv()
Return a changeable value of the pivot vector.
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition: DenseMatrix.h:71