3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MUTLIDOMAIN_COUPLING_JACOBIAN_PATTERN_HH
25#define DUMUX_MUTLIDOMAIN_COUPLING_JACOBIAN_PATTERN_HH
26
27#include <type_traits>
28#include <dune/common/indices.hh>
29#include <dune/istl/matrixindexset.hh>
31
32namespace Dumux {
33
39template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
40 typename std::enable_if_t<( (GridGeometryI::discMethod == DiscretizationMethods::cctpfa)
41 || (GridGeometryI::discMethod == DiscretizationMethods::ccmpfa) ), int> = 0>
42Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
43 Dune::index_constant<i> domainI,
44 const GridGeometryI& gridGeometryI,
45 Dune::index_constant<j> domainJ,
46 const GridGeometryJ& gridGeometryJ)
47{
48 const auto numDofsI = gridGeometryI.numDofs();
49 const auto numDofsJ = gridGeometryJ.numDofs();
50 Dune::MatrixIndexSet pattern;
51 pattern.resize(numDofsI, numDofsJ);
52
53 // matrix pattern for implicit Jacobians
54 if (isImplicit)
55 {
56 for (const auto& elementI : elements(gridGeometryI.gridView()))
57 {
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);
62 }
63 }
64
65 // matrix pattern for explicit Jacobians
66 // -> diagonal matrix, so coupling block is empty
67 // just return the empty pattern
68
69 return pattern;
70}
71
77template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
78 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::box), int> = 0>
79Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
80 Dune::index_constant<i> domainI,
81 const GridGeometryI& gridGeometryI,
82 Dune::index_constant<j> domainJ,
83 const GridGeometryJ& gridGeometryJ)
84{
85 const auto numDofsI = gridGeometryI.numDofs();
86 const auto numDofsJ = gridGeometryJ.numDofs();
87 Dune::MatrixIndexSet pattern;
88 pattern.resize(numDofsI, numDofsJ);
89
90 // matrix pattern for implicit Jacobians
91 if (isImplicit)
92 {
93 static constexpr int dim = std::decay_t<decltype(gridGeometryI.gridView())>::dimension;
94 for (const auto& elementI : elements(gridGeometryI.gridView()))
95 {
96 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, domainJ);
97 for (std::size_t vIdxLocal = 0; vIdxLocal < elementI.subEntities(dim); ++vIdxLocal)
98 {
99 const auto globalI = gridGeometryI.vertexMapper().subIndex(elementI, vIdxLocal, dim);
100 for (const auto globalJ : stencil)
101 pattern.add(globalI, globalJ);
102 }
103 }
104 }
105
106 // matrix pattern for explicit Jacobians
107 // -> diagonal matrix, so coupling block is empty
108 // just return the empty pattern
109
110 return pattern;
111}
112
118template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
119 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::staggered &&
120 GridGeometryI::isCellCenter()), int> = 0>
121Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
122 Dune::index_constant<i> domainI,
123 const GridGeometryI& gridGeometryI,
124 Dune::index_constant<j> domainJ,
125 const GridGeometryJ& gridGeometryJ)
126{
127 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
128
129 for (const auto& elementI : elements(gridGeometryI.gridView()))
130 {
131 const auto ccGlobalI = gridGeometryI.elementMapper().index(elementI);
132 for (auto&& faceGlobalJ : couplingManager.couplingStencil(domainI, elementI, domainJ))
133 pattern.add(ccGlobalI, faceGlobalJ);
134 }
135
136 return pattern;
137}
138
144template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
145 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::staggered &&
146 GridGeometryI::isFace()), int> = 0>
147Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
148 Dune::index_constant<i> domainI,
149 const GridGeometryI& gridGeometryI,
150 Dune::index_constant<j> domainJ,
151 const GridGeometryJ& gridGeometryJ)
152{
153 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
154
155 auto fvGeometry = localView(gridGeometryI);
156 for (const auto& elementI : elements(gridGeometryI.gridView()))
157 {
158 fvGeometry.bindElement(elementI);
159
160 // loop over sub control faces
161 for (auto&& scvf : scvfs(fvGeometry))
162 {
163 const auto faceGlobalI = scvf.dofIndex();
164 for (auto&& globalJ : couplingManager.couplingStencil(domainI, scvf, domainJ))
165 pattern.add(faceGlobalI, globalJ);
166 }
167 }
168
169 return pattern;
170}
171
177template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
178 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::fcstaggered), int> = 0>
179Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
180 Dune::index_constant<i> domainI,
181 const GridGeometryI& gridGeometryI,
182 Dune::index_constant<j> domainJ,
183 const GridGeometryJ& gridGeometryJ)
184{
185 Dune::MatrixIndexSet pattern(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
186
187 auto fvGeometry = localView(gridGeometryI);
188 for (const auto& elementI : elements(gridGeometryI.gridView()))
189 {
190 fvGeometry.bindElement(elementI);
191 for (const auto& scv : scvs(fvGeometry))
192 {
193 const auto globalI = scv.dofIndex();
194 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, scv, domainJ);
195 for (const auto globalJ : stencil)
196 {
197 assert(globalJ < gridGeometryJ.numDofs());
198 pattern.add(globalI, globalJ);
199
200 if (gridGeometryI.isPeriodic())
201 {
202 if (gridGeometryI.dofOnPeriodicBoundary(globalI))
203 {
204 const auto globalIP = gridGeometryI.periodicallyMappedDof(globalI);
205
206 if (globalI > globalIP)
207 pattern.add(globalIP, globalJ);
208 }
209 }
210 }
211 }
212 }
213
214 return pattern;
215}
216
222template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
223 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::fcdiamond), int> = 0>
224Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
225 Dune::index_constant<i> domainI,
226 const GridGeometryI& gridGeometryI,
227 Dune::index_constant<j> domainJ,
228 const GridGeometryJ& gridGeometryJ)
229{
230 Dune::MatrixIndexSet pattern;
231
232 // matrix pattern for implicit Jacobians
233 if (isImplicit)
234 {
235 pattern.resize(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
236 for (const auto& elementI : elements(gridGeometryI.gridView()))
237 {
238 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, domainJ);
239 for (std::size_t localFacetIndex = 0; localFacetIndex < elementI.subEntities(1); ++localFacetIndex)
240 {
241 const auto globalI = gridGeometryI.dofMapper().subIndex(elementI, localFacetIndex, 1);
242 for (const auto globalJ : stencil)
243 pattern.add(globalI, globalJ);
244 }
245 }
246 }
247
248 // matrix pattern for explicit Jacobians
249 // -> diagonal matrix, so coupling block is empty
250 // just return the empty pattern
251
252 return pattern;
253}
254
260template<bool isImplicit, class CouplingManager, class GridGeometryI, class GridGeometryJ, std::size_t i, std::size_t j,
261 typename std::enable_if_t<(GridGeometryI::discMethod == DiscretizationMethods::pq1bubble), int> = 0>
262Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager& couplingManager,
263 Dune::index_constant<i> domainI,
264 const GridGeometryI& gridGeometryI,
265 Dune::index_constant<j> domainJ,
266 const GridGeometryJ& gridGeometryJ)
267{
268 Dune::MatrixIndexSet pattern;
269
270 // matrix pattern for implicit Jacobians
271 if (isImplicit)
272 {
273 pattern.resize(gridGeometryI.numDofs(), gridGeometryJ.numDofs());
274 auto fvGeometry = localView(gridGeometryI);
275 for (const auto& elementI : elements(gridGeometryI.gridView()))
276 {
277 fvGeometry.bindElement(elementI);
278 const auto& stencil = couplingManager.couplingStencil(domainI, elementI, domainJ);
279 for (const auto& scv : scvs(fvGeometry))
280 {
281 for (const auto globalJ : stencil)
282 pattern.add(scv.dofIndex(), globalJ);
283
284 }
285 }
286 }
287
288 // matrix pattern for explicit Jacobians
289 // -> diagonal matrix, so coupling block is empty
290 // just return the empty pattern
291
292 return pattern;
293}
294
295} // end namespace Dumux
296
297#endif
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