24#ifndef DUMUX_MUTLIDOMAIN_COUPLING_JACOBIAN_PATTERN_HH
25#define DUMUX_MUTLIDOMAIN_COUPLING_JACOBIAN_PATTERN_HH
28#include <dune/common/indices.hh>
29#include <dune/istl/matrixindexset.hh>
39template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
43 Dune::index_constant<i> domainI,
44 const GridGeometryI& gridGeometryI,
45 Dune::index_constant<j> domainJ,
46 const GridGeometryJ& gridGeometryJ)
48 const auto numDofsI = gridGeometryI.numDofs();
49 const auto numDofsJ = gridGeometryJ.numDofs();
50 Dune::MatrixIndexSet pattern;
51 pattern.resize(numDofsI, numDofsJ);
56 for (
const auto& elementI : elements(gridGeometryI.gridView()))
58 const auto& stencil = couplingManager.
couplingStencil(domainI, elementI, domainJ);
59 const auto globalI = gridGeometryI.elementMapper().index(elementI);
60 for (
const auto globalJ : stencil)
61 pattern.add(globalI, globalJ);
77template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
80 Dune::index_constant<i> domainI,
81 const GridGeometryI& gridGeometryI,
82 Dune::index_constant<j> domainJ,
83 const GridGeometryJ& gridGeometryJ)
85 const auto numDofsI = gridGeometryI.numDofs();
86 const auto numDofsJ = gridGeometryJ.numDofs();
87 Dune::MatrixIndexSet pattern;
88 pattern.resize(numDofsI, numDofsJ);
93 static constexpr int dim = std::decay_t<
decltype(gridGeometryI.gridView())>::dimension;
94 for (
const auto& elementI : elements(gridGeometryI.gridView()))
96 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, domainJ);
97 for (std::size_t vIdxLocal = 0; vIdxLocal < elementI.subEntities(dim); ++vIdxLocal)
99 const auto globalI = gridGeometryI.vertexMapper().subIndex(elementI, vIdxLocal, dim);
100 for (
const auto globalJ : stencil)
101 pattern.add(globalI, globalJ);
118template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
120 GridGeometryI::isCellCenter()),
int> = 0>
122 Dune::index_constant<i> domainI,
123 const GridGeometryI& gridGeometryI,
124 Dune::index_constant<j> domainJ,
125 const GridGeometryJ& gridGeometryJ)
127 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
129 for (
const auto& elementI : elements(gridGeometryI.gridView()))
131 const auto ccGlobalI = gridGeometryI.elementMapper().index(elementI);
132 for (
auto&& faceGlobalJ : couplingManager.couplingStencil(domainI, elementI, domainJ))
133 pattern.add(ccGlobalI, faceGlobalJ);
144template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
146 GridGeometryI::isFace()),
int> = 0>
148 Dune::index_constant<i> domainI,
149 const GridGeometryI& gridGeometryI,
150 Dune::index_constant<j> domainJ,
151 const GridGeometryJ& gridGeometryJ)
153 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
155 for (
const auto& elementI : elements(gridGeometryI.gridView()))
157 auto fvGeometry =
localView(gridGeometryI);
158 fvGeometry.bindElement(elementI);
161 for (
auto&& scvf : scvfs(fvGeometry))
163 const auto faceGlobalI = scvf.dofIndex();
164 for (
auto&& globalJ : couplingManager.couplingStencil(domainI, scvf, domainJ))
165 pattern.add(faceGlobalI, globalJ);
The available discretization methods in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager &couplingManager, Dune::index_constant< i > domainI, const GridGeometryI &gridGeometryI, Dune::index_constant< j > domainJ, const GridGeometryJ &gridGeometryJ)
Helper function to generate coupling Jacobian pattern (off-diagonal blocks) for cell-centered schemes...
Definition: couplingjacobianpattern.hh:42
Definition: multidomain/couplingmanager.hh:46
const CouplingStencilType< i, j > & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &elementI, Dune::index_constant< j > domainJ) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: multidomain/couplingmanager.hh:83