12#ifndef DUMUX_MUTLIDOMAIN_COUPLING_JACOBIAN_PATTERN_HH
13#define DUMUX_MUTLIDOMAIN_COUPLING_JACOBIAN_PATTERN_HH
16#include <dune/common/indices.hh>
17#include <dune/istl/matrixindexset.hh>
28template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
32 Dune::index_constant<i> domainI,
33 const GridGeometryI& gridGeometryI,
34 Dune::index_constant<j> domainJ,
35 const GridGeometryJ& gridGeometryJ)
37 const auto numDofsI = gridGeometryI.numDofs();
38 const auto numDofsJ = gridGeometryJ.numDofs();
39 Dune::MatrixIndexSet pattern;
40 pattern.resize(numDofsI, numDofsJ);
45 for (
const auto& elementI : elements(gridGeometryI.gridView()))
47 const auto& stencil = couplingManager.
couplingStencil(domainI, elementI, domainJ);
48 const auto globalI = gridGeometryI.elementMapper().index(elementI);
49 for (
const auto globalJ : stencil)
50 pattern.add(globalI, globalJ);
66template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
68 GridGeometryI::isCellCenter()),
int> = 0>
70 Dune::index_constant<i> domainI,
71 const GridGeometryI& gridGeometryI,
72 Dune::index_constant<j> domainJ,
73 const GridGeometryJ& gridGeometryJ)
75 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
77 for (
const auto& elementI : elements(gridGeometryI.gridView()))
79 const auto ccGlobalI = gridGeometryI.elementMapper().index(elementI);
80 for (
auto&& faceGlobalJ : couplingManager.couplingStencil(domainI, elementI, domainJ))
81 pattern.add(ccGlobalI, faceGlobalJ);
92template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
94 GridGeometryI::isFace()),
int> = 0>
96 Dune::index_constant<i> domainI,
97 const GridGeometryI& gridGeometryI,
98 Dune::index_constant<j> domainJ,
99 const GridGeometryJ& gridGeometryJ)
101 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
103 auto fvGeometry =
localView(gridGeometryI);
104 for (
const auto& elementI : elements(gridGeometryI.gridView()))
106 fvGeometry.bindElement(elementI);
109 for (
auto&& scvf : scvfs(fvGeometry))
111 const auto faceGlobalI = scvf.dofIndex();
112 for (
auto&& globalJ : couplingManager.couplingStencil(domainI, scvf, domainJ))
113 pattern.add(faceGlobalI, globalJ);
125template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
128 Dune::index_constant<i> domainI,
129 const GridGeometryI& gridGeometryI,
130 Dune::index_constant<j> domainJ,
131 const GridGeometryJ& gridGeometryJ)
133 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
135 auto fvGeometry =
localView(gridGeometryI);
136 for (
const auto& elementI : elements(gridGeometryI.gridView()))
138 fvGeometry.bindElement(elementI);
139 for (
const auto& scv :
scvs(fvGeometry))
141 const auto globalI = scv.dofIndex();
142 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, scv, domainJ);
143 for (
const auto globalJ : stencil)
145 assert(globalJ < gridGeometryJ.numDofs());
146 pattern.add(globalI, globalJ);
148 if (gridGeometryI.isPeriodic())
150 if (gridGeometryI.dofOnPeriodicBoundary(globalI))
152 const auto globalIP = gridGeometryI.periodicallyMappedDof(globalI);
154 if (globalI > globalIP)
155 pattern.add(globalIP, globalJ);
170template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
171 typename std::enable_if_t<DiscretizationMethods::isCVFE<typename GridGeometryI::DiscretizationMethod>,
int> = 0>
173 Dune::index_constant<i> domainI,
174 const GridGeometryI& gridGeometryI,
175 Dune::index_constant<j> domainJ,
176 const GridGeometryJ& gridGeometryJ)
178 Dune::MatrixIndexSet pattern;
181 if constexpr (isImplicit)
183 pattern.resize(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
184 auto fvGeometry =
localView(gridGeometryI);
185 for (
const auto& elementI : elements(gridGeometryI.gridView()))
187 fvGeometry.bindElement(elementI);
188 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, domainJ);
189 for (
const auto& localDof :
localDofs(fvGeometry))
191 for (
const auto globalJ : stencil)
192 pattern.add(localDof.dofIndex(), globalJ);
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
const CouplingStencilType< i, j > & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &elementI, Dune::index_constant< j > domainJ) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: multidomain/couplingmanager.hh:103
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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:31
Class representing dofs on elements for control-volume finite element schemes.
The available discretization methods in Dumux.
constexpr CCMpfa ccmpfa
Definition: method.hh:146
constexpr CCTpfa cctpfa
Definition: method.hh:145
constexpr Staggered staggered
Definition: method.hh:149
constexpr FCStaggered fcstaggered
Definition: method.hh:151
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition: localdof.hh:79
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50