AutoDiffMatrix.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2014, 2015 SINTEF ICT, Applied Mathematics.
3 
4  This file is part of the Open Porous Media project (OPM).
5 
6  OPM is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  OPM is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with OPM. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #ifndef OPM_AUTODIFFMATRIX_HEADER_INCLUDED
21 #define OPM_AUTODIFFMATRIX_HEADER_INCLUDED
22 
23 #include <opm/common/utility/platform_dependent/disable_warnings.h>
24 
25 #include <Eigen/Eigen>
26 #include <Eigen/Sparse>
27 
28 #include <opm/common/utility/platform_dependent/reenable_warnings.h>
29 
30 #include <opm/common/ErrorMacros.hpp>
32 #include <vector>
33 
34 
35 namespace Opm
36 {
37 
44  {
45  public:
46  typedef std::vector<double> DiagRep;
47  typedef Eigen::SparseMatrix<double> SparseRep;
48 
49 
54  : type_(Zero),
55  rows_(0),
56  cols_(0),
57  diag_(),
58  sparse_()
59  {
60  }
61 
62 
66  AutoDiffMatrix(const int num_rows, const int num_cols)
67  : type_(Zero),
68  rows_(num_rows),
69  cols_(num_cols),
70  diag_(),
71  sparse_()
72  {
73  }
74 
75 
76 
77 
81  static AutoDiffMatrix createIdentity(const int num_rows_cols)
82  {
83  return AutoDiffMatrix(Identity, num_rows_cols, num_rows_cols);
84  }
85 
86 
87 
91  explicit AutoDiffMatrix(const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& d)
92  : type_(Diagonal),
93  rows_(d.rows()),
94  cols_(d.cols()),
95  diag_(d.diagonal().array().data(), d.diagonal().array().data() + d.rows()),
96  sparse_()
97  {
98  assert(rows_ == cols_);
99  }
100 
101 
102 
106  explicit AutoDiffMatrix(const Eigen::SparseMatrix<double>& s)
107  : type_(Sparse),
108  rows_(s.rows()),
109  cols_(s.cols()),
110  diag_(),
111  sparse_(s)
112  {
113  }
114 
115 
116 
117  AutoDiffMatrix(const AutoDiffMatrix& other) = default;
118  AutoDiffMatrix& operator=(const AutoDiffMatrix& other) = default;
119 
120 
121 
123  : type_(Zero),
124  rows_(0),
125  cols_(0),
126  diag_(),
127  sparse_()
128  {
129  swap(other);
130  }
131 
132 
133 
135  {
136  swap(other);
137  return *this;
138  }
139 
140 
141 
142  void swap(AutoDiffMatrix& other)
143  {
144  std::swap(type_, other.type_);
145  std::swap(rows_, other.rows_);
146  std::swap(cols_, other.cols_);
147  diag_.swap(other.diag_);
148  sparse_.swap(other.sparse_);
149  }
150 
151 
152 
160  {
161  assert(rows_ == rhs.rows_);
162  assert(cols_ == rhs.cols_);
163  switch (type_) {
164  case Zero:
165  return rhs;
166  case Identity:
167  switch (rhs.type_) {
168  case Zero:
169  return *this;
170  case Identity:
171  return addII(*this, rhs);
172  case Diagonal:
173  return addDI(rhs, *this);
174  case Sparse:
175  return addSI(rhs, *this);
176  default:
177  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << rhs.type_);
178  }
179  case Diagonal:
180  switch (rhs.type_) {
181  case Zero:
182  return *this;
183  case Identity:
184  return addDI(*this, rhs);
185  case Diagonal:
186  return addDD(*this, rhs);
187  case Sparse:
188  return addSD(rhs, *this);
189  default:
190  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << rhs.type_);
191  }
192  case Sparse:
193  switch (rhs.type_) {
194  case Zero:
195  return *this;
196  case Identity:
197  return addSI(*this, rhs);
198  case Diagonal:
199  return addSD(*this, rhs);
200  case Sparse:
201  return addSS(*this, rhs);
202  default:
203  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << rhs.type_);
204  }
205  default:
206  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << rhs.type_);
207  }
208  }
209 
210 
211 
212 
213 
214 
221  {
222  assert(cols_ == rhs.rows_);
223  switch (type_) {
224  case Zero:
225  return AutoDiffMatrix(rows_, rhs.cols_);
226  case Identity:
227  return rhs;
228  case Diagonal:
229  switch (rhs.type_) {
230  case Zero:
231  return AutoDiffMatrix(rows_, rhs.cols_);
232  case Identity:
233  return *this;
234  case Diagonal:
235  return mulDD(*this, rhs);
236  case Sparse:
237  return mulDS(*this, rhs);
238  default:
239  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << rhs.type_);
240  }
241  case Sparse:
242  switch (rhs.type_) {
243  case Zero:
244  return AutoDiffMatrix(rows_, rhs.cols_);
245  case Identity:
246  return *this;
247  case Diagonal:
248  return mulSD(*this, rhs);
249  case Sparse:
250  return mulSS(*this, rhs);
251  default:
252  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << rhs.type_);
253  }
254  default:
255  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << rhs.type_);
256  }
257  }
258 
259 
260 
261 
262 
263 
264 
266  {
267  *this = *this + rhs;
268  return *this;
269  }
270 
271 
272 
273 
274 
275 
277  {
278  *this = *this + (rhs * -1.0);
279  return *this;
280  }
281 
282 
283 
284 
285 
286 
292  AutoDiffMatrix operator*(const double rhs) const
293  {
294  switch (type_) {
295  case Zero:
296  return *this;
297  case Identity:
298  {
299  AutoDiffMatrix retval(*this);
300  retval.type_ = Diagonal;
301  retval.diag_.assign(rows_, rhs);
302  return retval;
303  }
304  case Diagonal:
305  {
306  AutoDiffMatrix retval(*this);
307  for (double& elem : retval.diag_) {
308  elem *= rhs;
309  }
310  return retval;
311  }
312  case Sparse:
313  {
314  AutoDiffMatrix retval(*this);
315  retval.sparse_ *= rhs;
316  return retval;
317  }
318  default:
319  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << type_);
320  }
321  }
322 
323 
324 
325 
326 
327 
333  AutoDiffMatrix operator/(const double rhs) const
334  {
335  switch (type_) {
336  case Zero:
337  return *this;
338  case Identity:
339  {
340  AutoDiffMatrix retval(*this);
341  retval.type_ = Diagonal;
342  retval.diag_.assign(rows_, 1.0/rhs);
343  return retval;
344  }
345  case Diagonal:
346  {
347  AutoDiffMatrix retval(*this);
348  for (double& elem : retval.diag_) {
349  elem /= rhs;
350  }
351  return retval;
352  }
353  case Sparse:
354  {
355  AutoDiffMatrix retval(*this);
356  retval.sparse_ /= rhs;
357  return retval;
358  }
359  default:
360  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << type_);
361  }
362  }
363 
364 
365 
366 
367 
368 
374  Eigen::VectorXd operator*(const Eigen::VectorXd& rhs) const
375  {
376  assert(cols_ == rhs.size());
377  switch (type_) {
378  case Zero:
379  return Eigen::VectorXd::Zero(rows_);
380  case Identity:
381  return rhs;
382  case Diagonal:
383  {
384  const Eigen::VectorXd d = Eigen::Map<const Eigen::VectorXd>(diag_.data(), rows_);
385  return d.cwiseProduct(rhs);
386  }
387  case Sparse:
388  return sparse_ * rhs;
389  default:
390  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << type_);
391  }
392  }
393 
394 
395 
396 
397 
398  // Add identity to identity
399  static AutoDiffMatrix addII(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
400  {
401  assert(lhs.type_ == Identity);
402  assert(rhs.type_ == Identity);
403  AutoDiffMatrix retval;
404  retval.type_ = Diagonal;
405  retval.rows_ = lhs.rows_;
406  retval.cols_ = rhs.cols_;
407  retval.diag_.assign(lhs.rows_, 2.0);
408  return retval;
409  }
410 
411  // Add diagonal to identity
412  static AutoDiffMatrix addDI(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
413  {
414  static_cast<void>(rhs); // Silence release-mode warning.
415  assert(lhs.type_ == Diagonal);
416  assert(rhs.type_ == Identity);
417  AutoDiffMatrix retval = lhs;
418  for (int r = 0; r < lhs.rows_; ++r) {
419  retval.diag_[r] += 1.0;
420  }
421  return retval;
422  }
423 
424  // Add diagonal to diagonal
425  static AutoDiffMatrix addDD(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
426  {
427  assert(lhs.type_ == Diagonal);
428  assert(rhs.type_ == Diagonal);
429  AutoDiffMatrix retval = lhs;
430  for (int r = 0; r < lhs.rows_; ++r) {
431  retval.diag_[r] += rhs.diag_[r];
432  }
433  return retval;
434  }
435 
436  // Add sparse to identity
437  static AutoDiffMatrix addSI(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
438  {
439  static_cast<void>(rhs); // Silence release-mode warning.
440  assert(lhs.type_ == Sparse);
441  assert(rhs.type_ == Identity);
442  AutoDiffMatrix retval = lhs;
443  retval.sparse_ += spdiag(Eigen::VectorXd::Ones(lhs.rows_));;
444  return retval;
445  }
446 
447  // Add sparse to diagonal
448  static AutoDiffMatrix addSD(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
449  {
450  assert(lhs.type_ == Sparse);
451  assert(rhs.type_ == Diagonal);
452  AutoDiffMatrix retval = lhs;
453  retval.sparse_ += spdiag(rhs.diag_);
454  return retval;
455  }
456 
457  // Add sparse to sparse
458  static AutoDiffMatrix addSS(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
459  {
460  assert(lhs.type_ == Sparse);
461  assert(rhs.type_ == Sparse);
462  AutoDiffMatrix retval = lhs;
463  retval.sparse_ += rhs.sparse_;
464  return retval;
465  }
466 
467 
468 
469 
470  // Multiply diagonal with diagonal
471  static AutoDiffMatrix mulDD(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
472  {
473  assert(lhs.type_ == Diagonal);
474  assert(rhs.type_ == Diagonal);
475  AutoDiffMatrix retval = lhs;
476  for (int r = 0; r < lhs.rows_; ++r) {
477  retval.diag_[r] *= rhs.diag_[r];
478  }
479  return retval;
480  }
481 
482  // Multiply diagonal with sparse
483  static AutoDiffMatrix mulDS(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
484  {
485  assert(lhs.type_ == Diagonal);
486  assert(rhs.type_ == Sparse);
487  AutoDiffMatrix retval;
488  retval.type_ = Sparse;
489  retval.rows_ = lhs.rows_;
490  retval.cols_ = rhs.cols_;
491  fastDiagSparseProduct(lhs.diag_, rhs.sparse_, retval.sparse_);
492  return retval;
493  }
494 
495  // Multiply sparse with diagonal
496  static AutoDiffMatrix mulSD(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
497  {
498  assert(lhs.type_ == Sparse);
499  assert(rhs.type_ == Diagonal);
500  AutoDiffMatrix retval;
501  retval.type_ = Sparse;
502  retval.rows_ = lhs.rows_;
503  retval.cols_ = rhs.cols_;
504  fastSparseDiagProduct(lhs.sparse_, rhs.diag_, retval.sparse_);
505  return retval;
506  }
507 
508  // Multiply sparse with sparse
509  static AutoDiffMatrix mulSS(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs)
510  {
511  assert(lhs.type_ == Sparse);
512  assert(rhs.type_ == Sparse);
513  AutoDiffMatrix retval;
514  retval.type_ = Sparse;
515  retval.rows_ = lhs.rows_;
516  retval.cols_ = rhs.cols_;
517  fastSparseProduct(lhs.sparse_, rhs.sparse_, retval.sparse_);
518  return retval;
519  }
520 
521 
522 
523 
529  template<class Scalar, int Options, class Index>
530  void toSparse(Eigen::SparseMatrix<Scalar, Options, Index>& s) const
531  {
532  switch (type_) {
533  case Zero:
534  s = Eigen::SparseMatrix<Scalar, Options, Index>(rows_, cols_);
535  return;
536  case Identity:
537  s = spdiag(Eigen::VectorXd::Ones(rows_));
538  return;
539  case Diagonal:
540  s = spdiag(diag_);
541  return;
542  case Sparse:
543  s = sparse_;
544  return;
545  }
546  }
547 
548 
549 
553  int rows() const
554  {
555  return rows_;
556  }
557 
558 
559 
563  int cols() const
564  {
565  return cols_;
566  }
567 
568 
569 
576  int nonZeros() const
577  {
578  switch (type_) {
579  case Zero:
580  return 0;
581  case Identity:
582  return rows_;
583  case Diagonal:
584  return rows_;
585  case Sparse:
586  return sparse_.nonZeros();
587  default:
588  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << type_);
589  }
590  }
591 
592 
593 
594 
598  double coeff(const int row, const int col) const
599  {
600  switch (type_) {
601  case Zero:
602  return 0.0;
603  case Identity:
604  return (row == col) ? 1.0 : 0.0;
605  case Diagonal:
606  return (row == col) ? diag_[row] : 0.0;
607  case Sparse:
608  return sparse_.coeff(row, col);
609  default:
610  OPM_THROW(std::logic_error, "Invalid AutoDiffMatrix type encountered: " << type_);
611  }
612  }
613 
614 
615 
616 
617 
623  const SparseRep& getSparse() const {
624  if (type_ != Sparse) {
632  SparseRep& mutable_sparse = const_cast<SparseRep&>(sparse_);
633  toSparse(mutable_sparse);
634  }
635  return sparse_;
636  }
637 
638 
639  private:
640  enum AudoDiffMatrixType { Zero, Identity, Diagonal, Sparse };
641 
642  AudoDiffMatrixType type_; //< Type of matrix
643  int rows_; //< Number of rows
644  int cols_; //< Number of columns
645  DiagRep diag_; //< Diagonal representation (only if type==Diagonal)
646  SparseRep sparse_; //< Sparse representation (only if type==Sparse)
647 
648 
649 
653  AutoDiffMatrix(AudoDiffMatrixType type, int rows, int cols,
654  DiagRep diag=DiagRep(), SparseRep sparse=SparseRep())
655  : type_(type),
656  rows_(rows),
657  cols_(cols),
658  diag_(diag),
659  sparse_(sparse)
660  {
661  }
662 
663 
664 
665 
666 
672  template <class V>
673  static inline
674  SparseRep
675  spdiag(const V& d)
676  {
677  const int n = d.size();
678  SparseRep mat(n, n);
679  mat.reserve(Eigen::ArrayXi::Ones(n, 1));
680  for (SparseRep::Index i = 0; i < n; ++i) {
681  if (d[i] != 0.0) {
682  mat.insert(i, i) = d[i];
683  }
684  }
685 
686  return mat;
687  }
688 
689  };
690 
691 
692 
696  inline void fastSparseProduct(const AutoDiffMatrix& lhs, const AutoDiffMatrix& rhs, AutoDiffMatrix& res)
697  {
698  res = lhs * rhs;
699  }
700 
701 
705  inline void fastSparseProduct(const Eigen::SparseMatrix<double>& lhs, const AutoDiffMatrix& rhs, AutoDiffMatrix& res)
706  {
707  res = AutoDiffMatrix(lhs) * rhs;
708  }
709 
710 
714  inline AutoDiffMatrix operator*(const Eigen::SparseMatrix<double>& lhs, const AutoDiffMatrix& rhs)
715  {
716  AutoDiffMatrix retval;
717  fastSparseProduct(lhs, rhs, retval);
718  return retval;
719  }
720 
721 } // namespace Opm
722 
723 
724 #endif // OPM_AUTODIFFMATRIX_HEADER_INCLUDED
AutoDiffMatrix operator*(const AutoDiffMatrix &rhs) const
Definition: AutoDiffMatrix.hpp:220
AutoDiff< Scalar > operator*(const AutoDiff< Scalar > &lhs, const AutoDiff< Scalar > &rhs)
Definition: AutoDiff.hpp:211
static AutoDiffMatrix mulDD(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:471
static AutoDiffMatrix createIdentity(const int num_rows_cols)
Definition: AutoDiffMatrix.hpp:81
Eigen::VectorXd operator*(const Eigen::VectorXd &rhs) const
Definition: AutoDiffMatrix.hpp:374
Definition: AdditionalObjectDeleter.hpp:22
static AutoDiffMatrix addDI(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:412
AutoDiffMatrix & operator+=(const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:265
int cols() const
Definition: AutoDiffMatrix.hpp:563
const SparseRep & getSparse() const
Definition: AutoDiffMatrix.hpp:623
void fastSparseDiagProduct(const Eigen::SparseMatrix< double > &lhs, const std::vector< double > &rhs, Eigen::SparseMatrix< double > &res)
Definition: fastSparseProduct.hpp:164
void toSparse(Eigen::SparseMatrix< Scalar, Options, Index > &s) const
Definition: AutoDiffMatrix.hpp:530
static AutoDiffMatrix mulSS(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:509
double coeff(const int row, const int col) const
Definition: AutoDiffMatrix.hpp:598
int nonZeros() const
Definition: AutoDiffMatrix.hpp:576
void fastSparseProduct(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs, AutoDiffMatrix &res)
Definition: AutoDiffMatrix.hpp:696
AutoDiffMatrix operator*(const double rhs) const
Definition: AutoDiffMatrix.hpp:292
Eigen::SparseMatrix< double > SparseRep
Definition: AutoDiffMatrix.hpp:47
AutoDiffMatrix & operator-=(const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:276
static AutoDiffMatrix addSD(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:448
void swap(AutoDiffMatrix &other)
Definition: AutoDiffMatrix.hpp:142
Definition: AutoDiffMatrix.hpp:43
AutoDiffMatrix operator+(const AutoDiffMatrix &rhs) const
Definition: AutoDiffMatrix.hpp:159
AutoDiffMatrix & operator=(const AutoDiffMatrix &other)=default
ADB::V V
Definition: BlackoilModelBase_impl.hpp:84
static AutoDiffMatrix mulSD(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:496
AutoDiffMatrix & operator=(AutoDiffMatrix &&other)
Definition: AutoDiffMatrix.hpp:134
std::vector< double > DiagRep
Definition: AutoDiffMatrix.hpp:46
AutoDiffMatrix(AutoDiffMatrix &&other)
Definition: AutoDiffMatrix.hpp:122
AutoDiffMatrix operator/(const double rhs) const
Definition: AutoDiffMatrix.hpp:333
static AutoDiffMatrix addSS(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:458
static AutoDiffMatrix addSI(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:437
AutoDiffMatrix()
Definition: AutoDiffMatrix.hpp:53
static AutoDiffMatrix mulDS(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:483
static AutoDiffMatrix addDD(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:425
static AutoDiffMatrix addII(const AutoDiffMatrix &lhs, const AutoDiffMatrix &rhs)
Definition: AutoDiffMatrix.hpp:399
AutoDiffMatrix(const Eigen::SparseMatrix< double > &s)
Definition: AutoDiffMatrix.hpp:106
AutoDiffMatrix(const int num_rows, const int num_cols)
Definition: AutoDiffMatrix.hpp:66
AutoDiffMatrix(const Eigen::DiagonalMatrix< double, Eigen::Dynamic > &d)
Definition: AutoDiffMatrix.hpp:91
void fastDiagSparseProduct(const std::vector< double > &lhs, const Eigen::SparseMatrix< double > &rhs, Eigen::SparseMatrix< double > &res)
Definition: fastSparseProduct.hpp:145
int rows() const
Definition: AutoDiffMatrix.hpp:553