36 using MatrixBlock =
typename Dune::FieldMatrix<Scalar, 1, 1>;
49 const auto numRows = getNumRows_(A);
52 auto M = BCRSMatrix(numRows, numRows, BCRSMatrix::random);
55 setOccupationPattern_(M, A);
69 static void setOccupationPattern_(BCRSMatrix& M,
const MultiTypeBlockMatrix& A)
72 const auto numRows = M.N();
73 Dune::MatrixIndexSet occupationPattern;
74 occupationPattern.resize(numRows, numRows);
77 auto addIndices = [&](
const auto& subMatrix,
const std::size_t startRow,
const std::size_t startCol)
80 static const Scalar eps = getParam<Scalar>(
"MatrixConverter.DeletePatternEntriesBelowAbsThreshold", -1.0);
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);
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)
98 std::size_t colIndex = 0;
99 forEach(A[i], [&](
const auto& subMatrix)
101 addIndices(subMatrix, rowIndex, colIndex);
103 using SubBlockType =
typename std::decay_t<
decltype(subMatrix)>::block_type;
105 colIndex += SubBlockType::cols * subMatrix.M();
108 if(colIndex == numRows)
109 rowIndex += SubBlockType::rows * subMatrix.N();
113 occupationPattern.exportIdx(M);
122 static void copyValues_(BCRSMatrix& M,
const MultiTypeBlockMatrix& A)
125 const auto numRows = M.N();
128 auto copyValues = [&](
const auto& subMatrix,
const std::size_t startRow,
const std::size_t startCol)
131 static const Scalar eps = getParam<Scalar>(
"MatrixConverter.DeletePatternEntriesBelowAbsThreshold", -1.0);
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];
145 using namespace Dune::Hybrid;
146 std::size_t rowIndex = 0;
147 forEach(std::make_index_sequence<MultiTypeBlockMatrix::N()>(), [&A, ©Values, &rowIndex, numRows](
const auto i)
149 std::size_t colIndex = 0;
150 forEach(A[i], [&](
const auto& subMatrix)
152 copyValues(subMatrix, rowIndex, colIndex);
154 using SubBlockType =
typename std::decay_t<
decltype(subMatrix)>::block_type;
156 colIndex += SubBlockType::cols * subMatrix.M();
159 if(colIndex == numRows)
160 rowIndex += SubBlockType::rows * subMatrix.N();
170 static std::size_t getNumRows_(
const MultiTypeBlockMatrix& A)
173 std::size_t numRows = 0;
174 Dune::Hybrid::forEach(Dune::Hybrid::elementAt(A, Dune::Indices::_0), [&numRows](
const auto& subMatrixInFirstRow)
177 const auto numEq = std::decay_t<
decltype(subMatrixInFirstRow)>::block_type::cols;
178 numRows += numEq * subMatrixInFirstRow.M();