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 auto fvGeometry =
localView(gridGeometryI);
156 for (
const auto& elementI : elements(gridGeometryI.gridView()))
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);
177template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
180 Dune::index_constant<i> domainI,
181 const GridGeometryI& gridGeometryI,
182 Dune::index_constant<j> domainJ,
183 const GridGeometryJ& gridGeometryJ)
185 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
187 auto fvGeometry =
localView(gridGeometryI);
188 for (
const auto& elementI : elements(gridGeometryI.gridView()))
190 fvGeometry.bindElement(elementI);
191 for (
const auto& scv : scvs(fvGeometry))
193 const auto globalI = scv.dofIndex();
194 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, scv, domainJ);
195 for (
const auto globalJ : stencil)
197 assert(globalJ < gridGeometryJ.numDofs());
198 pattern.add(globalI, globalJ);
200 if (gridGeometryI.isPeriodic())
202 if (gridGeometryI.dofOnPeriodicBoundary(globalI))
204 const auto globalIP = gridGeometryI.periodicallyMappedDof(globalI);
206 if (globalI > globalIP)
207 pattern.add(globalIP, globalJ);
222template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
225 Dune::index_constant<i> domainI,
226 const GridGeometryI& gridGeometryI,
227 Dune::index_constant<j> domainJ,
228 const GridGeometryJ& gridGeometryJ)
230 Dune::MatrixIndexSet pattern;
235 pattern.resize(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
236 for (
const auto& elementI : elements(gridGeometryI.gridView()))
238 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, domainJ);
239 for (std::size_t localFacetIndex = 0; localFacetIndex < elementI.subEntities(1); ++localFacetIndex)
241 const auto globalI = gridGeometryI.dofMapper().subIndex(elementI, localFacetIndex, 1);
242 for (
const auto globalJ : stencil)
243 pattern.add(globalI, globalJ);
260template<
bool isImplicit,
class CouplingManager,
class GridGeometryI,
class GridGeometryJ, std::size_t i, std::size_t j,
263 Dune::index_constant<i> domainI,
264 const GridGeometryI& gridGeometryI,
265 Dune::index_constant<j> domainJ,
266 const GridGeometryJ& gridGeometryJ)
268 Dune::MatrixIndexSet pattern;
273 pattern.resize(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
274 auto fvGeometry =
localView(gridGeometryI);
275 for (
const auto& elementI : elements(gridGeometryI.gridView()))
277 fvGeometry.bindElement(elementI);
278 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, domainJ);
279 for (
const auto& scv : scvs(fvGeometry))
281 for (
const auto globalJ : stencil)
282 pattern.add(scv.dofIndex(), 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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr CCMpfa ccmpfa
Definition: method.hh:135
constexpr FCDiamond fcdiamond
Definition: method.hh:141
constexpr CCTpfa cctpfa
Definition: method.hh:134
constexpr Box box
Definition: method.hh:136
constexpr Staggered staggered
Definition: method.hh:138
constexpr PQ1Bubble pq1bubble
Definition: method.hh:137
constexpr FCStaggered fcstaggered
Definition: method.hh:140
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
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:126