version 3.11-dev
couplingjacobianpattern.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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_MUTLIDOMAIN_COUPLING_JACOBIAN_PATTERN_HH
13#define DUMUX_MUTLIDOMAIN_COUPLING_JACOBIAN_PATTERN_HH
14
15#include <type_traits>
16#include <dune/common/indices.hh>
17#include <dune/istl/matrixindexset.hh>
20
21namespace Dumux {
22
28template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
29 typename std::enable_if_t<( (GridGeometryI::discMethod == DiscretizationMethods::cctpfa)
30 || (GridGeometryI::discMethod == DiscretizationMethods::ccmpfa) ), int> = 0>
31Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
32 Dune::index_constant<i> domainI,
33 const GridGeometryI& gridGeometryI,
34 Dune::index_constant<j> domainJ,
35 const GridGeometryJ& gridGeometryJ)
36{
37 const auto numDofsI = gridGeometryI.numDofs();
38 const auto numDofsJ = gridGeometryJ.numDofs();
39 Dune::MatrixIndexSet pattern;
40 pattern.resize(numDofsI, numDofsJ);
41
42 // matrix pattern for implicit Jacobians
43 if (isImplicit)
44 {
45 for (const auto& elementI : elements(gridGeometryI.gridView()))
46 {
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);
51 }
52 }
53
54 // matrix pattern for explicit Jacobians
55 // -> diagonal matrix, so coupling block is empty
56 // just return the empty pattern
57
58 return pattern;
59}
60
66template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
67 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::staggered &&
68 GridGeometryI::isCellCenter()), int> = 0>
69Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
70 Dune::index_constant<i> domainI,
71 const GridGeometryI& gridGeometryI,
72 Dune::index_constant<j> domainJ,
73 const GridGeometryJ& gridGeometryJ)
74{
75 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
76
77 for (const auto& elementI : elements(gridGeometryI.gridView()))
78 {
79 const auto ccGlobalI = gridGeometryI.elementMapper().index(elementI);
80 for (auto&& faceGlobalJ : couplingManager.couplingStencil(domainI, elementI, domainJ))
81 pattern.add(ccGlobalI, faceGlobalJ);
82 }
83
84 return pattern;
85}
86
92template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
93 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::staggered &&
94 GridGeometryI::isFace()), int> = 0>
95Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
96 Dune::index_constant<i> domainI,
97 const GridGeometryI& gridGeometryI,
98 Dune::index_constant<j> domainJ,
99 const GridGeometryJ& gridGeometryJ)
100{
101 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
102
103 auto fvGeometry = localView(gridGeometryI);
104 for (const auto& elementI : elements(gridGeometryI.gridView()))
105 {
106 fvGeometry.bindElement(elementI);
107
108 // loop over sub control faces
109 for (auto&& scvf : scvfs(fvGeometry))
110 {
111 const auto faceGlobalI = scvf.dofIndex();
112 for (auto&& globalJ : couplingManager.couplingStencil(domainI, scvf, domainJ))
113 pattern.add(faceGlobalI, globalJ);
114 }
115 }
116
117 return pattern;
118}
119
125template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
126 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::fcstaggered), int> = 0>
127Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
128 Dune::index_constant<i> domainI,
129 const GridGeometryI& gridGeometryI,
130 Dune::index_constant<j> domainJ,
131 const GridGeometryJ& gridGeometryJ)
132{
133 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
134
135 auto fvGeometry = localView(gridGeometryI);
136 for (const auto& elementI : elements(gridGeometryI.gridView()))
137 {
138 fvGeometry.bindElement(elementI);
139 for (const auto& scv : scvs(fvGeometry))
140 {
141 const auto globalI = scv.dofIndex();
142 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, scv, domainJ);
143 for (const auto globalJ : stencil)
144 {
145 assert(globalJ < gridGeometryJ.numDofs());
146 pattern.add(globalI, globalJ);
147
148 if (gridGeometryI.isPeriodic())
149 {
150 if (gridGeometryI.dofOnPeriodicBoundary(globalI))
151 {
152 const auto globalIP = gridGeometryI.periodicallyMappedDof(globalI);
153
154 if (globalI > globalIP)
155 pattern.add(globalIP, globalJ);
156 }
157 }
158 }
159 }
160 }
161
162 return pattern;
163}
164
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>
172Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
173 Dune::index_constant<i> domainI,
174 const GridGeometryI& gridGeometryI,
175 Dune::index_constant<j> domainJ,
176 const GridGeometryJ& gridGeometryJ)
177{
178 Dune::MatrixIndexSet pattern;
179
180 // matrix pattern for implicit Jacobians
181 if constexpr (isImplicit)
182 {
183 pattern.resize(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
184 auto fvGeometry = localView(gridGeometryI);
185 for (const auto& elementI : elements(gridGeometryI.gridView()))
186 {
187 fvGeometry.bindElement(elementI);
188 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, domainJ);
189 for (const auto& localDof : localDofs(fvGeometry))
190 {
191 for (const auto globalJ : stencil)
192 pattern.add(localDof.dofIndex(), globalJ);
193
194 }
195 }
196 }
197
198 // matrix pattern for explicit Jacobians
199 // -> diagonal matrix, so coupling block is empty
200 // just return the empty pattern
201
202 return pattern;
203}
204
205} // end namespace Dumux
206
207#endif
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
Definition: adapt.hh:17
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