DuneMatrix.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2014 SINTEF ICT, Applied Mathematics.
3  Copyright 2014 IRIS AS
4 
5  This file is part of the Open Porous Media project (OPM).
6 
7  OPM is free software: you can redistribute it and/or modify
8  it under the terms of the GNU General Public License as published by
9  the Free Software Foundation, either version 3 of the License, or
10  (at your option) any later version.
11 
12  OPM is distributed in the hope that it will be useful,
13  but WITHOUT ANY WARRANTY; without even the implied warranty of
14  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  GNU General Public License for more details.
16 
17  You should have received a copy of the GNU General Public License
18  along with OPM. If not, see <http://www.gnu.org/licenses/>.
19 */
20 
21 #ifndef OPM_DUNEMATRIX_HEADER_INCLUDED
22 #define OPM_DUNEMATRIX_HEADER_INCLUDED
23 
24 #ifdef DUNE_BCRSMATRIX_HH
25 #error This header must be included before any bcrsmatrix.hh is included (directly or indirectly)
26 #endif
27 
28 #include <opm/common/utility/platform_dependent/disable_warnings.h>
29 
30 #include <Eigen/Eigen>
31 #include <Eigen/Sparse>
32 
33 #include <dune/common/fmatrix.hh>
34 #include <dune/common/version.hh>
35 
36 #if DUNE_VERSION_NEWER(DUNE_ISTL,2,4)
37 #include <dune/istl/bcrsmatrix.hh>
38 #else
39 // Include matrix header with hackery to make it possible to inherit.
40 #define private protected
41 #include <dune/istl/bcrsmatrix.hh>
42 #undef private
43 #endif
44 
45 #include <opm/common/utility/platform_dependent/reenable_warnings.h>
46 
47 namespace Opm
48 {
49 
50  class DuneMatrix : public Dune::BCRSMatrix< Dune::FieldMatrix<double, 1, 1> >
51  {
52  public:
53  DuneMatrix(const int rows, const int cols, const int* ia, const int* ja, const double* sa)
54  {
55  // create BCRSMatrix from given CSR storage
56  init( rows, cols, ia, ja, sa );
57  }
58 
60  DuneMatrix( const Eigen::SparseMatrix<double, Eigen::RowMajor>& matrix )
61  {
62  // Create ISTL matrix.
63  const int rows = matrix.rows();
64  const int cols = matrix.cols();
65  const int* ia = matrix.outerIndexPtr();
66  const int* ja = matrix.innerIndexPtr();
67  const double* sa = matrix.valuePtr();
68  // create BCRSMatrix from Eigen matrix
69  init( rows, cols, ia, ja, sa );
70  }
71 
72  protected:
73  void init(const int rows, const int cols, const int* ia, const int* ja, const double* sa)
74  {
75  typedef Dune::BCRSMatrix< Dune::FieldMatrix<double, 1, 1> > Super;
76  typedef Super::block_type block_type;
77  this->build_mode = Super::unknown;
78  this->ready = Super::built;
79  this->n = rows;
80  this->m = cols;
81  this->nnz = ia[rows];
82 
83 #if DUNE_VERSION_NEWER(DUNE_ISTL, 2, 3)
84  this->allocationSize = this->nnz;
85  this->avg = 0;
86  this->overflowsize = -1.0;
87 #endif
88 
89  // make sure to use the allocators of this matrix
90  // because the same allocators are used to deallocate the data
91  this->a = this->allocator_.allocate(this->nnz);
92  static_assert(sizeof(block_type) == sizeof(double), "This constructor requires a block type that is the same as a double.");
93  std::copy(sa, sa + this->nnz, reinterpret_cast<double*>(this->a));
94  this->j.reset(this->sizeAllocator_.allocate(this->nnz));
95  std::copy(ja, ja +this-> nnz, this->j.get());
96  this->r = rowAllocator_.allocate(rows);
97  for (int row = 0; row < rows; ++row) {
98  this->r[row].set(ia[row+1] - ia[row], this->a + ia[row], this->j.get() + ia[row]);
99  }
100  }
101  };
102 
103 } // namespace Opm
104 
105 #endif // OPM_DUNEMATRIX_HEADER_INCLUDED
Definition: AdditionalObjectDeleter.hpp:22
void init(const int rows, const int cols, const int *ia, const int *ja, const double *sa)
Definition: DuneMatrix.hpp:73
DuneMatrix(const int rows, const int cols, const int *ia, const int *ja, const double *sa)
Definition: DuneMatrix.hpp:53
Definition: DuneMatrix.hpp:50
DuneMatrix(const Eigen::SparseMatrix< double, Eigen::RowMajor > &matrix)
create an ISTL BCRSMatrix from a Eigen::SparseMatrix
Definition: DuneMatrix.hpp:60