12#ifndef DUMUX_JACOBIAN_PATTERN_HH
13#define DUMUX_JACOBIAN_PATTERN_HH
16#include <dune/istl/matrixindexset.hh>
25template<
bool isImplicit,
class GridGeometry,
30 const auto numDofs = gridGeometry.numDofs();
31 Dune::MatrixIndexSet pattern;
32 pattern.resize(numDofs, numDofs);
37 for (
unsigned int globalI = 0; globalI < numDofs; ++globalI)
39 pattern.add(globalI, globalI);
40 for (
const auto& dataJ : gridGeometry.connectivityMap()[globalI])
41 pattern.add(dataJ.globalJ, globalI);
48 for (
unsigned int globalI = 0; globalI < numDofs; ++globalI)
49 pattern.add(globalI, globalI);
59template<
bool isImplicit,
class GridGeometry,
64 const auto numDofs = gridGeometry.numDofs();
65 Dune::MatrixIndexSet pattern(numDofs, numDofs);
67 const auto& connectivityMap = gridGeometry.connectivityMap();
69 auto fvGeometry =
localView(gridGeometry);
71 for (
const auto& element : elements(gridGeometry.gridView()))
73 if(gridGeometry.isCellCenter())
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);
84 static constexpr auto faceIdx = GridGeometry::faceIdx();
85 fvGeometry.bindElement(element);
88 for (
auto&& scvf : scvfs(fvGeometry))
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);
105template<
class FEBasis>
108 const auto numDofs = feBasis.size();
110 Dune::MatrixIndexSet pattern;
111 pattern.resize(numDofs, numDofs);
115 for (
const auto& element : elements(feBasis.gridView()))
119 const auto& finiteElement =
localView.tree().finiteElement();
120 const auto numLocalDofs = finiteElement.localBasis().size();
121 for (std::size_t i = 0; i < numLocalDofs; i++)
124 for (std::size_t j = 0; j < numLocalDofs; j++)
127 pattern.add(dofIdxI, dofIdxJ);
141template<
bool isImplicit,
class GridGeometry,
150template<
bool isImplicit,
class GridGeometry,
155 const auto numDofs = gridGeometry.numDofs();
156 Dune::MatrixIndexSet pattern(numDofs, numDofs);
158 const auto& connectivityMap = gridGeometry.connectivityMap();
159 auto fvGeometry =
localView(gridGeometry);
162 for (
const auto& element : elements(gridGeometry.gridView()))
164 fvGeometry.bind(element);
165 for (
const auto& scv : scvs(fvGeometry))
167 const auto globalI = scv.dofIndex();
168 pattern.add(globalI, globalI);
170 for (
const auto& scvIdxJ : connectivityMap[scv.index()])
172 const auto globalJ = fvGeometry.scv(scvIdxJ).dofIndex();
173 pattern.add(globalI, globalJ);
175 if (gridGeometry.isPeriodic())
177 if (gridGeometry.dofOnPeriodicBoundary(globalI) && globalI != globalJ)
179 const auto globalIP = gridGeometry.periodicallyMappedDof(globalI);
180 pattern.add(globalIP, globalI);
181 pattern.add(globalI, globalIP);
183 if (globalI > globalIP)
184 pattern.add(globalIP, globalJ);
198template<
bool isImplicit,
class GridGeometry,
199 typename std::enable_if_t<DiscretizationMethods::isCVFE<typename GridGeometry::DiscretizationMethod>,
int> = 0>
203 const auto numDofs = gridGeometry.numDofs();
204 Dune::MatrixIndexSet pattern(numDofs, numDofs);
207 if constexpr (isImplicit)
209 auto fvGeometry =
localView(gridGeometry);
210 for (
const auto& element : elements(gridGeometry.gridView()))
212 fvGeometry.bindElement(element);
213 for (
const auto& scv : scvs(fvGeometry))
215 const auto globalI = scv.dofIndex();
216 for (
const auto& scvJ : scvs(fvGeometry))
218 const auto globalJ = scvJ.dofIndex();
219 pattern.add(globalI, globalJ);
221 if (gridGeometry.isPeriodic())
223 if (gridGeometry.dofOnPeriodicBoundary(globalI) && globalI != globalJ)
225 const auto globalIP = gridGeometry.periodicallyMappedDof(globalI);
226 pattern.add(globalIP, globalI);
227 pattern.add(globalI, globalIP);
229 if (globalI > globalIP)
230 pattern.add(globalIP, globalJ);
241 for (
unsigned int globalI = 0; globalI < numDofs; ++globalI)
242 pattern.add(globalI, globalI);
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