21 #ifndef OPM_DUNEMATRIX_HEADER_INCLUDED
22 #define OPM_DUNEMATRIX_HEADER_INCLUDED
24 #ifdef DUNE_BCRSMATRIX_HH
25 #error This header must be included before any bcrsmatrix.hh is included (directly or indirectly)
28 #include <opm/common/utility/platform_dependent/disable_warnings.h>
30 #include <Eigen/Eigen>
31 #include <Eigen/Sparse>
33 #include <dune/common/fmatrix.hh>
34 #include <dune/common/version.hh>
36 #include <dune/istl/bcrsmatrix.hh>
38 #include <opm/common/utility/platform_dependent/reenable_warnings.h>
43 class DuneMatrix :
public Dune::BCRSMatrix< Dune::FieldMatrix<double, 1, 1> >
46 DuneMatrix(
const int rows,
const int cols,
const int* ia,
const int* ja,
const double* sa)
49 init( rows, cols, ia, ja, sa );
53 DuneMatrix(
const Eigen::SparseMatrix<double, Eigen::RowMajor>& matrix )
56 const int rows = matrix.rows();
57 const int cols = matrix.cols();
58 const int* ia = matrix.outerIndexPtr();
59 const int* ja = matrix.innerIndexPtr();
60 const double* sa = matrix.valuePtr();
62 init( rows, cols, ia, ja, sa );
66 void init(
const int rows,
const int cols,
const int* ia,
const int* ja,
const double* sa)
68 typedef Dune::BCRSMatrix< Dune::FieldMatrix<double, 1, 1> > Super;
69 typedef Super::block_type block_type;
70 this->build_mode = Super::unknown;
71 this->ready = Super::built;
75 typedef Super::size_type size_type ;
76 #if DUNE_VERSION_NEWER_REV(DUNE_ISTL, 2, 4, 1)
77 size_type& nnz = this->nnz_;
78 std::shared_ptr<size_type>& j = this->j_;
80 size_type& nnz = this->nnz;
81 std::shared_ptr<size_type>& j = this->j;
86 #if DUNE_VERSION_NEWER(DUNE_ISTL, 2, 5)
87 this->allocationSize_ = nnz;
89 this->allocationSize = nnz;
92 this->overflowsize = -1.0;
96 this->a = this->allocator_.allocate(nnz);
97 static_assert(
sizeof(block_type) ==
sizeof(
double),
"This constructor requires a block type that is the same as a double.");
98 std::copy(sa, sa + nnz, reinterpret_cast<double*>(this->a));
99 j.reset(this->sizeAllocator_.allocate(nnz));
100 std::copy(ja, ja +nnz, j.get());
101 this->r = rowAllocator_.allocate(rows);
102 for (
int row = 0; row < rows; ++row) {
103 this->r[row].set(ia[row+1] - ia[row], this->a + ia[row], j.get() + ia[row]);
110 #endif // OPM_DUNEMATRIX_HEADER_INCLUDED
DuneMatrix(const Eigen::SparseMatrix< double, Eigen::RowMajor > &matrix)
create an ISTL BCRSMatrix from a Eigen::SparseMatrix
Definition: DuneMatrix.hpp:53
Definition: DuneMatrix.hpp:43