istlsparsematrixadapter.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
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 2 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 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef EWOMS_ISTL_SPARSE_MATRIX_ADAPTER_HH
28#define EWOMS_ISTL_SPARSE_MATRIX_ADAPTER_HH
29
30#include <dune/istl/bcrsmatrix.hh>
31#include <dune/common/fmatrix.hh>
32#include <dune/common/version.hh>
33
34namespace Opm {
35namespace Linear {
36
41template <class MatrixBlockType, class AllocatorType=std::allocator<MatrixBlockType> >
43{
44public:
46 using IstlMatrix = Dune::BCRSMatrix<MatrixBlockType, AllocatorType>;
47
49 using MatrixBlock = typename IstlMatrix::block_type;
50 static_assert(std::is_same<MatrixBlock, MatrixBlockType>::value,
51 "IstlMatrix::block_type and MatrixBlockType must be identical");
52
54 using Scalar = typename MatrixBlock::field_type;
55
59 IstlSparseMatrixAdapter(const size_t rows, const size_t columns)
60 : rows_(rows)
61 , columns_(columns)
62 , istlMatrix_()
63 {}
64
68 template <class Simulator>
70 : IstlSparseMatrixAdapter(simulator.model().numTotalDof(), simulator.model().numTotalDof())
71 {}
72
76 template <class Set>
77 void reserve(const std::vector<Set>& sparsityPattern)
78 {
79 // allocate raw matrix
80 istlMatrix_.reset(new IstlMatrix(rows_, columns_, IstlMatrix::random));
81
82 // make sure sparsityPattern is consistent with number of rows
83 assert(rows_ == sparsityPattern.size());
84
85 // allocate space for the rows of the matrix
86 for (size_t dofIdx = 0; dofIdx < rows_; ++ dofIdx)
87 istlMatrix_->setrowsize(dofIdx, sparsityPattern[dofIdx].size());
88
89 istlMatrix_->endrowsizes();
90
91 // fill the rows with indices. each degree of freedom talks to
92 // all of its neighbors. (it also talks to itself since
93 // degrees of freedom are sometimes quite egocentric.)
94 for (size_t dofIdx = 0; dofIdx < rows_; ++ dofIdx) {
95 auto nIt = sparsityPattern[dofIdx].begin();
96 auto nEndIt = sparsityPattern[dofIdx].end();
97 for (; nIt != nEndIt; ++nIt)
98 istlMatrix_->addindex(dofIdx, *nIt);
99 }
100 istlMatrix_->endindices();
101 }
102
107 { return *istlMatrix_; }
108 const IstlMatrix& istlMatrix() const
109 { return *istlMatrix_; }
110
114 size_t rows() const
115 { return rows_; }
116
120 size_t cols() const
121 { return columns_; }
122
126 void clear()
127 { (*istlMatrix_) = Scalar(0.0); }
128
135 void clearRow(const size_t row, const Scalar diag = 1.0)
136 {
137 MatrixBlock diagBlock(Scalar(0));
138 for (int i = 0; i < diagBlock.rows; ++i)
139 diagBlock[i][i] = diag;
140
141 auto& matRow = (*istlMatrix_)[row];
142 auto colIt = matRow.begin();
143 const auto& colEndIt = matRow.end();
144 for (; colIt != colEndIt; ++colIt) {
145 if (colIt.index() == row)
146 *colIt = diagBlock;
147 else
148 *colIt = Scalar(0.0);
149 }
150 }
151
155 void block(const size_t rowIdx, const size_t colIdx, MatrixBlock& value) const
156 { value = (*istlMatrix_)[rowIdx][colIdx]; }
157
158 MatrixBlockType* blockAddress(const size_t rowIdx, const size_t colIdx) const
159 { return &(*istlMatrix_)[rowIdx][colIdx]; }
160
161
165 void setBlock(const size_t rowIdx, const size_t colIdx, const MatrixBlock& value)
166 { (*istlMatrix_)[rowIdx][colIdx] = value; }
167
171 void addToBlock(const size_t rowIdx, const size_t colIdx, const MatrixBlock& value)
172 { (*istlMatrix_)[rowIdx][colIdx] += value; }
173
179 void commit()
180 { }
181
189 void finalize()
190 { }
191
192protected:
193 size_t rows_;
194 size_t columns_;
195
196 std::unique_ptr<IstlMatrix> istlMatrix_;
197};
198
199}} // namespace Linear, Opm
200
201#endif
A sparse matrix interface backend for BCRSMatrix from dune-istl.
Definition: istlsparsematrixadapter.hh:43
size_t columns_
Definition: istlsparsematrixadapter.hh:194
void setBlock(const size_t rowIdx, const size_t colIdx, const MatrixBlock &value)
Set matrix block to given block.
Definition: istlsparsematrixadapter.hh:165
void clear()
Set all matrix entries to zero.
Definition: istlsparsematrixadapter.hh:126
Dune::BCRSMatrix< MatrixBlockType, AllocatorType > IstlMatrix
Implementation of matrix.
Definition: istlsparsematrixadapter.hh:46
void finalize()
Finish modifying the matrix, i.e., convert the data structure from one tailored for linearization to ...
Definition: istlsparsematrixadapter.hh:189
std::unique_ptr< IstlMatrix > istlMatrix_
Definition: istlsparsematrixadapter.hh:196
IstlSparseMatrixAdapter(const size_t rows, const size_t columns)
Constructor creating an empty matrix.
Definition: istlsparsematrixadapter.hh:59
void reserve(const std::vector< Set > &sparsityPattern)
Allocate matrix structure give a sparsity pattern.
Definition: istlsparsematrixadapter.hh:77
size_t rows_
Definition: istlsparsematrixadapter.hh:193
size_t cols() const
Return number of columns of the matrix.
Definition: istlsparsematrixadapter.hh:120
size_t rows() const
Return number of rows of the matrix.
Definition: istlsparsematrixadapter.hh:114
const IstlMatrix & istlMatrix() const
Definition: istlsparsematrixadapter.hh:108
typename IstlMatrix::block_type MatrixBlock
block type forming the matrix entries
Definition: istlsparsematrixadapter.hh:49
void clearRow(const size_t row, const Scalar diag=1.0)
Set given row to zero except for the main-diagonal entry (if it exists).
Definition: istlsparsematrixadapter.hh:135
void commit()
Commit matrix from local caches into matrix native structure.
Definition: istlsparsematrixadapter.hh:179
void addToBlock(const size_t rowIdx, const size_t colIdx, const MatrixBlock &value)
Add block to matrix block.
Definition: istlsparsematrixadapter.hh:171
IstlMatrix & istlMatrix()
Return constant reference to matrix implementation.
Definition: istlsparsematrixadapter.hh:106
void block(const size_t rowIdx, const size_t colIdx, MatrixBlock &value) const
Fill given block with entries stored in the matrix.
Definition: istlsparsematrixadapter.hh:155
typename MatrixBlock::field_type Scalar
type of scalar
Definition: istlsparsematrixadapter.hh:54
IstlSparseMatrixAdapter(const Simulator &simulator)
Constructor taking simulator and creating an empty matrix .
Definition: istlsparsematrixadapter.hh:69
MatrixBlockType * blockAddress(const size_t rowIdx, const size_t colIdx) const
Definition: istlsparsematrixadapter.hh:158
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:102
Definition: blackoilboundaryratevector.hh:37