24#ifndef DUMUX_JACOBIAN_PATTERN_HH
25#define DUMUX_JACOBIAN_PATTERN_HH
28#include <dune/istl/matrixindexset.hh>
37template<
bool isImplicit,
class GridGeometry,
41 const auto numDofs = gridGeometry.numDofs();
42 Dune::MatrixIndexSet pattern;
43 pattern.resize(numDofs, numDofs);
48 static constexpr int dim = std::decay_t<
decltype(gridGeometry.gridView())>::dimension;
49 for (
const auto& element : elements(gridGeometry.gridView()))
51 for (
unsigned int vIdx = 0; vIdx < element.subEntities(dim); ++vIdx)
53 const auto globalI = gridGeometry.vertexMapper().subIndex(element, vIdx, dim);
54 for (
unsigned int vIdx2 = 0; vIdx2 < element.subEntities(dim); ++vIdx2)
56 const auto globalJ = gridGeometry.vertexMapper().subIndex(element, vIdx2, dim);
57 pattern.add(globalI, globalJ);
59 if (gridGeometry.dofOnPeriodicBoundary(globalI) && globalI != globalJ)
61 const auto globalIP = gridGeometry.periodicallyMappedDof(globalI);
62 pattern.add(globalIP, globalI);
63 pattern.add(globalI, globalIP);
64 if (globalI > globalIP)
65 pattern.add(globalIP, globalJ);
75 for (
unsigned int globalI = 0; globalI < numDofs; ++globalI)
76 pattern.add(globalI, globalI);
86template<
bool isImplicit,
class GridGeometry,
91 const auto numDofs = gridGeometry.numDofs();
92 Dune::MatrixIndexSet pattern;
93 pattern.resize(numDofs, numDofs);
98 for (
unsigned int globalI = 0; globalI < numDofs; ++globalI)
100 pattern.add(globalI, globalI);
101 for (
const auto& dataJ : gridGeometry.connectivityMap()[globalI])
102 pattern.add(dataJ.globalJ, globalI);
109 for (
unsigned int globalI = 0; globalI < numDofs; ++globalI)
110 pattern.add(globalI, globalI);
120template<
bool isImplicit,
class GridGeometry,
125 const auto numDofs = gridGeometry.numDofs();
126 Dune::MatrixIndexSet pattern(numDofs, numDofs);
128 const auto& connectivityMap = gridGeometry.connectivityMap();
131 for (
const auto& element : elements(gridGeometry.gridView()))
133 if(gridGeometry.isCellCenter())
136 static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
137 const auto ccGlobalI = gridGeometry.elementMapper().index(element);
138 pattern.add(ccGlobalI, ccGlobalI);
139 for (
auto&& ccGlobalJ : connectivityMap(cellCenterIdx, cellCenterIdx, ccGlobalI))
140 pattern.add(ccGlobalI, ccGlobalJ);
144 static constexpr auto faceIdx = GridGeometry::faceIdx();
145 auto fvGeometry =
localView(gridGeometry);
146 fvGeometry.bindElement(element);
149 for (
auto&& scvf : scvfs(fvGeometry))
151 const auto faceGlobalI = scvf.dofIndex();
152 pattern.add(faceGlobalI, faceGlobalI);
153 for (
auto&& faceGlobalJ : connectivityMap(faceIdx, faceIdx, scvf.index()))
154 pattern.add(faceGlobalI, faceGlobalJ);
166template<
class FEBasis>
169 const auto numDofs = feBasis.size();
171 Dune::MatrixIndexSet pattern;
172 pattern.resize(numDofs, numDofs);
175 for (
const auto& element : elements(feBasis.gridView()))
180 const auto& finiteElement =
localView.tree().finiteElement();
181 const auto numLocalDofs = finiteElement.localBasis().size();
182 for (
size_t i = 0; i < numLocalDofs; i++)
185 for (
size_t j = 0; j < numLocalDofs; j++)
188 pattern.add(dofIdxI, dofIdxJ);
202template<
bool isImplicit,
class GridGeometry,
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 getJacobianPattern(const GridGeometry &gridGeometry)
Helper function to generate Jacobian pattern for the box method.
Definition: jacobianpattern.hh:39
Dune::MatrixIndexSet getFEJacobianPattern(const FEBasis &feBasis)
Helper function to generate Jacobian pattern for finite element scheme.
Definition: jacobianpattern.hh:167
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29