3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MATRIX_CONVERTER
25#define DUMUX_MATRIX_CONVERTER
26
27#include <cmath>
28#include <dune/common/indices.hh>
29#include <dune/common/hybridutilities.hh>
30#include <dune/istl/bvector.hh>
31#include <dune/istl/bcrsmatrix.hh>
32#include <dune/istl/matrixindexset.hh>
33#include <dune/istl/multitypeblockvector.hh>
34#include <dune/istl/multitypeblockmatrix.hh>
35
37
38namespace Dumux {
39
46template <class MultiTypeBlockMatrix, class Scalar=double>
48{
49 using MatrixBlock = typename Dune::FieldMatrix<Scalar, 1, 1>;
50 using BCRSMatrix = typename Dune::BCRSMatrix<MatrixBlock>;
51
52public:
53
59 static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
60 {
61 // get the size for the converted matrix
62 const auto numRows = getNumRows_(A);
63
64 // create an empty BCRS matrix with 1x1 blocks
65 auto M = BCRSMatrix(numRows, numRows, BCRSMatrix::random);
66
67 // set the occupation pattern and copy the values
68 setOccupationPattern_(M, A);
69 copyValues_(M, A);
70
71 return M;
72 }
73
74private:
75
82 static void setOccupationPattern_(BCRSMatrix& M, const MultiTypeBlockMatrix& A)
83 {
84 // prepare the occupation pattern
85 const auto numRows = M.N();
86 Dune::MatrixIndexSet occupationPattern;
87 occupationPattern.resize(numRows, numRows);
88
89 // lambda function to fill the occupation pattern
90 auto addIndices = [&occupationPattern](const auto& subMatrix, const std::size_t startRow, const std::size_t startCol)
91 {
92 using std::abs;
93 static const Scalar eps = getParam<Scalar>("MatrixConverter.DeletePatternEntriesBelowAbsThreshold", -1.0);
94
95 using BlockType = typename std::decay_t<decltype(subMatrix)>::block_type;
96 const auto blockSizeI = BlockType::rows;
97 const auto blockSizeJ = BlockType::cols;
98 for(auto row = subMatrix.begin(); row != subMatrix.end(); ++row)
99 for(auto col = row->begin(); col != row->end(); ++col)
100 for(std::size_t i = 0; i < blockSizeI; ++i)
101 for(std::size_t j = 0; j < blockSizeJ; ++j)
102 if(abs(subMatrix[row.index()][col.index()][i][j]) > eps)
103 occupationPattern.add(startRow + row.index()*blockSizeI + i, startCol + col.index()*blockSizeJ + j);
104 };
105
106 // fill the pattern
107 std::size_t rowIndex = 0;
108 Dune::Hybrid::forEach(A, [&addIndices, &rowIndex, numRows](const auto& rowOfMultiTypeMatrix)
109 {
110 std::size_t colIndex = 0;
111 Dune::Hybrid::forEach(rowOfMultiTypeMatrix, [&addIndices, &colIndex, &rowIndex, numRows](const auto& subMatrix)
112 {
113 addIndices(subMatrix, rowIndex, colIndex);
114
115 using SubBlockType = typename std::decay_t<decltype(subMatrix)>::block_type;
116
117 colIndex += SubBlockType::cols * subMatrix.M();
118
119 // if we have arrived at the right side of the matrix, increase the row index
120 if(colIndex == numRows)
121 rowIndex += SubBlockType::rows * subMatrix.N();
122 });
123 });
124
125 occupationPattern.exportIdx(M);
126 }
127
134 static void copyValues_(BCRSMatrix& M, const MultiTypeBlockMatrix& A)
135 {
136 // get number of rows
137 const auto numRows = M.N();
138
139 // lambda function to copy the values
140 auto copyValues = [&M](const auto& subMatrix, const std::size_t startRow, const std::size_t startCol)
141 {
142 using std::abs;
143 static const Scalar eps = getParam<Scalar>("MatrixConverter.DeletePatternEntriesBelowAbsThreshold", -1.0);
144
145 using BlockType = typename std::decay_t<decltype(subMatrix)>::block_type;
146 const auto blockSizeI = BlockType::rows;
147 const auto blockSizeJ = BlockType::cols;
148 for (auto row = subMatrix.begin(); row != subMatrix.end(); ++row)
149 for (auto col = row->begin(); col != row->end(); ++col)
150 for (std::size_t i = 0; i < blockSizeI; ++i)
151 for (std::size_t j = 0; j < blockSizeJ; ++j)
152 if(abs(subMatrix[row.index()][col.index()][i][j]) > eps)
153 M[startRow + row.index()*blockSizeI + i][startCol + col.index()*blockSizeJ + j] = subMatrix[row.index()][col.index()][i][j];
154
155 };
156
157 std::size_t rowIndex = 0;
158 Dune::Hybrid::forEach(A, [&copyValues, &rowIndex, numRows](const auto& rowOfMultiTypeMatrix)
159 {
160 std::size_t colIndex = 0;
161 Dune::Hybrid::forEach(rowOfMultiTypeMatrix, [&copyValues, &colIndex, &rowIndex, numRows](const auto& subMatrix)
162 {
163 copyValues(subMatrix, rowIndex, colIndex);
164
165 using SubBlockType = typename std::decay_t<decltype(subMatrix)>::block_type;
166
167 colIndex += SubBlockType::cols * subMatrix.M();
168
169 // if we have arrived at the right side of the matrix, increase the row index
170 if(colIndex == numRows)
171 rowIndex += SubBlockType::rows * subMatrix.N();
172 });
173 });
174 }
175
181 static std::size_t getNumRows_(const MultiTypeBlockMatrix& A)
182 {
183 // iterate over the first row of the multitype blockmatrix
184 std::size_t numRows = 0;
185 Dune::Hybrid::forEach(Dune::Hybrid::elementAt(A, Dune::Indices::_0), [&numRows](const auto& subMatrixInFirstRow)
186 {
187 // the number of cols of the individual submatrice's block equals the respective number of equations.
188 const auto numEq = std::decay_t<decltype(subMatrixInFirstRow)>::block_type::cols;
189 numRows += numEq * subMatrixInFirstRow.M();
190 });
191
192 return numRows;
193 }
194
195};
196
201template<class MultiTypeBlockVector, class Scalar=double>
203{
204 using VectorBlock = typename Dune::FieldVector<Scalar, 1>;
205 using BlockVector = typename Dune::BlockVector<VectorBlock>;
206
207public:
208
214 static auto multiTypeToBlockVector(const MultiTypeBlockVector& b)
215 {
216 const auto size = getSize_(b);
217
218 BlockVector bTmp;
219 bTmp.resize(size);
220
221 std::size_t startIndex = 0;
222 Dune::Hybrid::forEach(b, [&bTmp, &startIndex](const auto& subVector)
223 {
224 const auto numEq = std::decay_t<decltype(subVector)>::block_type::size();
225
226 for(std::size_t i = 0; i < subVector.size(); ++i)
227 for(std::size_t j = 0; j < numEq; ++j)
228 bTmp[startIndex + i*numEq + j] = subVector[i][j];
229
230 startIndex += numEq*subVector.size();
231 });
232
233 return bTmp;
234 }
235
242 static void retrieveValues(MultiTypeBlockVector& x, const BlockVector& y)
243 {
244 std::size_t startIndex = 0;
245 Dune::Hybrid::forEach(x, [&y, &startIndex](auto& subVector)
246 {
247 const auto numEq = std::decay_t<decltype(subVector)>::block_type::size();
248
249 for(std::size_t i = 0; i < subVector.size(); ++i)
250 for(std::size_t j = 0; j < numEq; ++j)
251 subVector[i][j] = y[startIndex + i*numEq + j];
252
253 startIndex += numEq*subVector.size();
254 });
255 }
256
257private:
258
264 static std::size_t getSize_(const MultiTypeBlockVector& b)
265 {
266 std::size_t size = 0;
267 Dune::Hybrid::forEach(b, [&size](const auto& subVector)
268 {
269 // the size of the individual vector blocks equals the respective number of equations.
270 const auto numEq = std::decay_t<decltype(subVector)>::block_type::size();
271 size += numEq * subVector.size();
272 });
273
274 return size;
275 }
276};
277
278} // end namespace Dumux
279
280#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
A helper classe that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix TODO: allow ...
Definition: matrixconverter.hh:48
static auto multiTypeToBCRSMatrix(const MultiTypeBlockMatrix &A)
Converts the matrix to a type the IterativeSolverBackend can handle.
Definition: matrixconverter.hh:59
A helper classe that converts a Dune::MultiTypeBlockVector into a plain Dune::BlockVector and transfe...
Definition: matrixconverter.hh:203
static void retrieveValues(MultiTypeBlockVector &x, const BlockVector &y)
Copys the entries of a Dune::BlockVector to a Dune::MultiTypeBlockVector.
Definition: matrixconverter.hh:242
static auto multiTypeToBlockVector(const MultiTypeBlockVector &b)
Converts a Dune::MultiTypeBlockVector to a plain 1x1 Dune::BlockVector.
Definition: matrixconverter.hh:214