version 3.10-dev
matrixconverter.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_MATRIX_CONVERTER
13#define DUMUX_MATRIX_CONVERTER
14
15#include <cmath>
16#include <utility>
17#include <dune/common/indices.hh>
18#include <dune/common/hybridutilities.hh>
19#include <dune/istl/bvector.hh>
20#include <dune/istl/bcrsmatrix.hh>
21#include <dune/istl/matrixindexset.hh>
22
24
25namespace Dumux {
26
33template <class MultiTypeBlockMatrix, class Scalar=double>
35{
36 using MatrixBlock = typename Dune::FieldMatrix<Scalar, 1, 1>;
37 using BCRSMatrix = typename Dune::BCRSMatrix<MatrixBlock>;
38
39public:
40
46 static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
47 {
48 // get the size for the converted matrix
49 const auto numRows = getNumRows_(A);
50
51 // create an empty BCRS matrix with 1x1 blocks
52 auto M = BCRSMatrix(numRows, numRows, BCRSMatrix::random);
53
54 // set the occupation pattern and copy the values
55 setOccupationPattern_(M, A);
56 copyValues_(M, A);
57
58 return M;
59 }
60
61private:
62
69 static void setOccupationPattern_(BCRSMatrix& M, const MultiTypeBlockMatrix& A)
70 {
71 // prepare the occupation pattern
72 const auto numRows = M.N();
73 Dune::MatrixIndexSet occupationPattern;
74 occupationPattern.resize(numRows, numRows);
75
76 // lambda function to fill the occupation pattern
77 auto addIndices = [&](const auto& subMatrix, const std::size_t startRow, const std::size_t startCol)
78 {
79 using std::abs;
80 static const Scalar eps = getParam<Scalar>("MatrixConverter.DeletePatternEntriesBelowAbsThreshold", -1.0);
81
82 using BlockType = typename std::decay_t<decltype(subMatrix)>::block_type;
83 const auto blockSizeI = BlockType::rows;
84 const auto blockSizeJ = BlockType::cols;
85 for(auto row = subMatrix.begin(); row != subMatrix.end(); ++row)
86 for(auto col = row->begin(); col != row->end(); ++col)
87 for(std::size_t i = 0; i < blockSizeI; ++i)
88 for(std::size_t j = 0; j < blockSizeJ; ++j)
89 if(abs(subMatrix[row.index()][col.index()][i][j]) > eps)
90 occupationPattern.add(startRow + row.index()*blockSizeI + i, startCol + col.index()*blockSizeJ + j);
91 };
92
93 // fill the pattern
94 using namespace Dune::Hybrid;
95 std::size_t rowIndex = 0;
96 forEach(std::make_index_sequence<MultiTypeBlockMatrix::N()>(), [&A, &addIndices, &rowIndex, numRows](const auto i)
97 {
98 std::size_t colIndex = 0;
99 forEach(A[i], [&](const auto& subMatrix)
100 {
101 addIndices(subMatrix, rowIndex, colIndex);
102
103 using SubBlockType = typename std::decay_t<decltype(subMatrix)>::block_type;
104
105 colIndex += SubBlockType::cols * subMatrix.M();
106
107 // if we have arrived at the right side of the matrix, increase the row index
108 if(colIndex == numRows)
109 rowIndex += SubBlockType::rows * subMatrix.N();
110 });
111 });
112
113 occupationPattern.exportIdx(M);
114 }
115
122 static void copyValues_(BCRSMatrix& M, const MultiTypeBlockMatrix& A)
123 {
124 // get number of rows
125 const auto numRows = M.N();
126
127 // lambda function to copy the values
128 auto copyValues = [&](const auto& subMatrix, const std::size_t startRow, const std::size_t startCol)
129 {
130 using std::abs;
131 static const Scalar eps = getParam<Scalar>("MatrixConverter.DeletePatternEntriesBelowAbsThreshold", -1.0);
132
133 using BlockType = typename std::decay_t<decltype(subMatrix)>::block_type;
134 const auto blockSizeI = BlockType::rows;
135 const auto blockSizeJ = BlockType::cols;
136 for (auto row = subMatrix.begin(); row != subMatrix.end(); ++row)
137 for (auto col = row->begin(); col != row->end(); ++col)
138 for (std::size_t i = 0; i < blockSizeI; ++i)
139 for (std::size_t j = 0; j < blockSizeJ; ++j)
140 if(abs(subMatrix[row.index()][col.index()][i][j]) > eps)
141 M[startRow + row.index()*blockSizeI + i][startCol + col.index()*blockSizeJ + j] = subMatrix[row.index()][col.index()][i][j];
142
143 };
144
145 using namespace Dune::Hybrid;
146 std::size_t rowIndex = 0;
147 forEach(std::make_index_sequence<MultiTypeBlockMatrix::N()>(), [&A, &copyValues, &rowIndex, numRows](const auto i)
148 {
149 std::size_t colIndex = 0;
150 forEach(A[i], [&](const auto& subMatrix)
151 {
152 copyValues(subMatrix, rowIndex, colIndex);
153
154 using SubBlockType = typename std::decay_t<decltype(subMatrix)>::block_type;
155
156 colIndex += SubBlockType::cols * subMatrix.M();
157
158 // if we have arrived at the right side of the matrix, increase the row index
159 if(colIndex == numRows)
160 rowIndex += SubBlockType::rows * subMatrix.N();
161 });
162 });
163 }
164
170 static std::size_t getNumRows_(const MultiTypeBlockMatrix& A)
171 {
172 // iterate over the first row of the multitype blockmatrix
173 std::size_t numRows = 0;
174 Dune::Hybrid::forEach(Dune::Hybrid::elementAt(A, Dune::Indices::_0), [&numRows](const auto& subMatrixInFirstRow)
175 {
176 // the number of cols of the individual submatrice's block equals the respective number of equations.
177 const auto numEq = std::decay_t<decltype(subMatrixInFirstRow)>::block_type::cols;
178 numRows += numEq * subMatrixInFirstRow.M();
179 });
180
181 return numRows;
182 }
183
184};
185
190template<class MultiTypeBlockVector, class Scalar=double>
192{
193 using VectorBlock = typename Dune::FieldVector<Scalar, 1>;
194 using BlockVector = typename Dune::BlockVector<VectorBlock>;
195
196public:
197
203 static auto multiTypeToBlockVector(const MultiTypeBlockVector& b)
204 {
205 BlockVector bTmp;
206 bTmp.resize(b.dim());
207
208 std::size_t startIndex = 0;
209 Dune::Hybrid::forEach(b, [&](const auto& subVector)
210 {
211 const auto numEq = std::decay_t<decltype(subVector)>::block_type::size();
212
213 for(std::size_t i = 0; i < subVector.size(); ++i)
214 for(std::size_t j = 0; j < numEq; ++j)
215 bTmp[startIndex + i*numEq + j] = subVector[i][j];
216
217 startIndex += numEq*subVector.size();
218 });
219
220 return bTmp;
221 }
222
229 static void retrieveValues(MultiTypeBlockVector& x, const BlockVector& y)
230 {
231 std::size_t startIndex = 0;
232 Dune::Hybrid::forEach(x, [&](auto& subVector)
233 {
234 const auto numEq = std::decay_t<decltype(subVector)>::block_type::size();
235
236 for(std::size_t i = 0; i < subVector.size(); ++i)
237 for(std::size_t j = 0; j < numEq; ++j)
238 subVector[i][j] = y[startIndex + i*numEq + j];
239
240 startIndex += numEq*subVector.size();
241 });
242 }
243};
244
245} // end namespace Dumux
246
247#endif
A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix TODO: allow b...
Definition: matrixconverter.hh:35
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition: matrixconverter.hh:46
A helper class that converts a Dune::MultiTypeBlockVector into a plain Dune::BlockVector and transfer...
Definition: matrixconverter.hh:192
static void retrieveValues(MultiTypeBlockVector &x, const BlockVector &y)
Copies the entries of a Dune::BlockVector to a Dune::MultiTypeBlockVector.
Definition: matrixconverter.hh:229
static auto multiTypeToBlockVector(const MultiTypeBlockVector &b)
Converts a Dune::MultiTypeBlockVector to a plain 1x1 Dune::BlockVector.
Definition: matrixconverter.hh:203
Definition: matrix.hh:20
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.