version 3.8
jacobianpattern.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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_JACOBIAN_PATTERN_HH
13#define DUMUX_JACOBIAN_PATTERN_HH
14
15#include <type_traits>
16#include <dune/istl/matrixindexset.hh>
18
19namespace Dumux {
20
25template<bool isImplicit, class GridGeometry,
26 typename std::enable_if_t<( (GridGeometry::discMethod == DiscretizationMethods::cctpfa)
27 || (GridGeometry::discMethod == DiscretizationMethods::ccmpfa) ), int> = 0>
28Dune::MatrixIndexSet getJacobianPattern(const GridGeometry& gridGeometry)
29{
30 const auto numDofs = gridGeometry.numDofs();
31 Dune::MatrixIndexSet pattern;
32 pattern.resize(numDofs, numDofs);
33
34 // matrix pattern for implicit Jacobians
35 if (isImplicit)
36 {
37 for (unsigned int globalI = 0; globalI < numDofs; ++globalI)
38 {
39 pattern.add(globalI, globalI);
40 for (const auto& dataJ : gridGeometry.connectivityMap()[globalI])
41 pattern.add(dataJ.globalJ, globalI);
42 }
43 }
44
45 // matrix pattern for explicit Jacobians -> diagonal matrix
46 else
47 {
48 for (unsigned int globalI = 0; globalI < numDofs; ++globalI)
49 pattern.add(globalI, globalI);
50 }
51
52 return pattern;
53}
54
59template<bool isImplicit, class GridGeometry,
60 typename std::enable_if_t<( (GridGeometry::discMethod == DiscretizationMethods::staggered) ), int> = 0>
61auto getJacobianPattern(const GridGeometry& gridGeometry)
62{
63 // resize the jacobian and the residual
64 const auto numDofs = gridGeometry.numDofs();
65 Dune::MatrixIndexSet pattern(numDofs, numDofs);
66
67 const auto& connectivityMap = gridGeometry.connectivityMap();
68
69 auto fvGeometry = localView(gridGeometry);
70 // evaluate the actual pattern
71 for (const auto& element : elements(gridGeometry.gridView()))
72 {
73 if(gridGeometry.isCellCenter())
74 {
75 // the global index of the element at hand
76 static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
77 const auto ccGlobalI = gridGeometry.elementMapper().index(element);
78 pattern.add(ccGlobalI, ccGlobalI);
79 for (auto&& ccGlobalJ : connectivityMap(cellCenterIdx, cellCenterIdx, ccGlobalI))
80 pattern.add(ccGlobalI, ccGlobalJ);
81 }
82 else
83 {
84 static constexpr auto faceIdx = GridGeometry::faceIdx();
85 fvGeometry.bindElement(element);
86
87 // loop over sub control faces
88 for (auto&& scvf : scvfs(fvGeometry))
89 {
90 const auto faceGlobalI = scvf.dofIndex();
91 pattern.add(faceGlobalI, faceGlobalI);
92 for (auto&& faceGlobalJ : connectivityMap(faceIdx, faceIdx, scvf.index()))
93 pattern.add(faceGlobalI, faceGlobalJ);
94 }
95 }
96 }
97
98 return pattern;
99}
100
105template<class FEBasis>
106Dune::MatrixIndexSet getFEJacobianPattern(const FEBasis& feBasis)
107{
108 const auto numDofs = feBasis.size();
109
110 Dune::MatrixIndexSet pattern;
111 pattern.resize(numDofs, numDofs);
112
113 auto localView = feBasis.localView();
114 // matrix pattern for implicit Jacobians
115 for (const auto& element : elements(feBasis.gridView()))
116 {
117 localView.bind(element);
118
119 const auto& finiteElement = localView.tree().finiteElement();
120 const auto numLocalDofs = finiteElement.localBasis().size();
121 for (std::size_t i = 0; i < numLocalDofs; i++)
122 {
123 const auto dofIdxI = localView.index(i);
124 for (std::size_t j = 0; j < numLocalDofs; j++)
125 {
126 const auto dofIdxJ = localView.index(j);
127 pattern.add(dofIdxI, dofIdxJ);
128 }
129 }
130 }
131
132 return pattern;
133}
134
141template<bool isImplicit, class GridGeometry,
142 typename std::enable_if_t<(GridGeometry::discMethod == DiscretizationMethods::fem), int> = 0>
143Dune::MatrixIndexSet getJacobianPattern(const GridGeometry& gridGeometry)
144{ return getFEJacobianPattern(gridGeometry.feBasis()); }
145
150template<bool isImplicit, class GridGeometry,
151 typename std::enable_if_t<(GridGeometry::discMethod == DiscretizationMethods::fcstaggered), int> = 0>
152Dune::MatrixIndexSet getJacobianPattern(const GridGeometry& gridGeometry)
153{
154 // resize the jacobian and the residual
155 const auto numDofs = gridGeometry.numDofs();
156 Dune::MatrixIndexSet pattern(numDofs, numDofs);
157
158 const auto& connectivityMap = gridGeometry.connectivityMap();
159 auto fvGeometry = localView(gridGeometry);
160
161 // set the pattern
162 for (const auto& element : elements(gridGeometry.gridView()))
163 {
164 fvGeometry.bind(element);
165 for (const auto& scv : scvs(fvGeometry))
166 {
167 const auto globalI = scv.dofIndex();
168 pattern.add(globalI, globalI);
169
170 for (const auto& scvIdxJ : connectivityMap[scv.index()])
171 {
172 const auto globalJ = fvGeometry.scv(scvIdxJ).dofIndex();
173 pattern.add(globalI, globalJ);
174
175 if (gridGeometry.isPeriodic())
176 {
177 if (gridGeometry.dofOnPeriodicBoundary(globalI) && globalI != globalJ)
178 {
179 const auto globalIP = gridGeometry.periodicallyMappedDof(globalI);
180 pattern.add(globalIP, globalI);
181 pattern.add(globalI, globalIP);
182
183 if (globalI > globalIP)
184 pattern.add(globalIP, globalJ);
185 }
186 }
187 }
188 }
189 }
190
191 return pattern;
192}
193
198template<bool isImplicit, class GridGeometry,
199 typename std::enable_if_t<DiscretizationMethods::isCVFE<typename GridGeometry::DiscretizationMethod>, int> = 0>
200Dune::MatrixIndexSet getJacobianPattern(const GridGeometry& gridGeometry)
201{
202 // resize the jacobian and the residual
203 const auto numDofs = gridGeometry.numDofs();
204 Dune::MatrixIndexSet pattern(numDofs, numDofs);
205
206 // matrix pattern for implicit Jacobians
207 if constexpr (isImplicit)
208 {
209 auto fvGeometry = localView(gridGeometry);
210 for (const auto& element : elements(gridGeometry.gridView()))
211 {
212 fvGeometry.bindElement(element);
213 for (const auto& scv : scvs(fvGeometry))
214 {
215 const auto globalI = scv.dofIndex();
216 for (const auto& scvJ : scvs(fvGeometry))
217 {
218 const auto globalJ = scvJ.dofIndex();
219 pattern.add(globalI, globalJ);
220
221 if (gridGeometry.isPeriodic())
222 {
223 if (gridGeometry.dofOnPeriodicBoundary(globalI) && globalI != globalJ)
224 {
225 const auto globalIP = gridGeometry.periodicallyMappedDof(globalI);
226 pattern.add(globalIP, globalI);
227 pattern.add(globalI, globalIP);
228
229 if (globalI > globalIP)
230 pattern.add(globalIP, globalJ);
231 }
232 }
233 }
234 }
235 }
236 }
237
238 // matrix pattern for explicit Jacobians -> diagonal matrix
239 else
240 {
241 for (unsigned int globalI = 0; globalI < numDofs; ++globalI)
242 pattern.add(globalI, globalI);
243 }
244
245 return pattern;
246}
247
248
249} // namespace Dumux
250
251#endif
Dune::MatrixIndexSet getJacobianPattern(const GridGeometry &gridGeometry)
Helper function to generate Jacobian pattern for cell-centered methods.
Definition: jacobianpattern.hh:28
Dune::MatrixIndexSet getFEJacobianPattern(const FEBasis &feBasis)
Helper function to generate Jacobian pattern for finite element scheme.
Definition: jacobianpattern.hh:106
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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 FEM fem
Definition: method.hh:150
constexpr FCStaggered fcstaggered
Definition: method.hh:151
Definition: adapt.hh:17