12#ifndef DUMUX_DISCRETIZATION_FLUXSTENCIL_HH
13#define DUMUX_DISCRETIZATION_FLUXSTENCIL_HH
17#include <dune/common/reservedvector.hh>
32template<
class FVElementGeometry,
class DiscretizationMethod =
typename FVElementGeometry::Gr
idGeometry::DiscretizationMethod>
40template<
class FVElementGeometry>
41class FluxStencil<FVElementGeometry, DiscretizationMethods::CCTpfa>
43 using GridGeometry =
typename FVElementGeometry::GridGeometry;
44 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
45 using GridView =
typename GridGeometry::GridView;
46 using Element =
typename GridView::template Codim<0>::Entity;
54 using Stencil =
typename SubControlVolumeFace::Traits::GridIndexStorage;
58 const FVElementGeometry& fvGeometry,
59 const SubControlVolumeFace& scvf)
62 return Stencil({scvf.insideScvIdx()});
63 else if (scvf.numOutsideScvs() > 1)
65 Stencil stencil({scvf.insideScvIdx()});
66 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
67 stencil.push_back(scvf.outsideScvIdx(i));
71 return Stencil({scvf.insideScvIdx(), scvf.outsideScvIdx()});
80template<
class FVElementGeometry>
81class FluxStencil<FVElementGeometry, DiscretizationMethods::CCMpfa>
83 using GridGeometry =
typename FVElementGeometry::GridGeometry;
84 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
85 using GridView =
typename GridGeometry::GridView;
86 using Element =
typename GridView::template Codim<0>::Entity;
90 using NodalIndexSet =
typename GridGeometry::GridIVIndexSets::DualGridIndexSet::NodalIndexSet;
97 using Stencil =
typename NodalIndexSet::NodalGridStencilType;
101 const FVElementGeometry& fvGeometry,
102 const SubControlVolumeFace& scvf)
104 const auto& gridGeometry = fvGeometry.gridGeometry();
107 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
108 return gridGeometry.gridInteractionVolumeIndexSets().secondaryIndexSet(scvf).gridScvIndices();
110 return gridGeometry.gridInteractionVolumeIndexSets().primaryIndexSet(scvf).gridScvIndices();
std::vector< GridIndexType > ScvfStencilIForJ
We don't know yet how many faces couple to a neighboring element.
Definition: fluxstencil.hh:94
static const Stencil & stencil(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Returns the indices of the elements required for flux calculation on an scvf.
Definition: fluxstencil.hh:100
typename NodalIndexSet::NodalGridStencilType Stencil
The flux stencil type.
Definition: fluxstencil.hh:97
Dune::ReservedVector< GridIndexType, 1 > ScvfStencilIForJ
Each cell I couples to a cell J always only via one face.
Definition: fluxstencil.hh:51
typename SubControlVolumeFace::Traits::GridIndexStorage Stencil
The flux stencil type.
Definition: fluxstencil.hh:54
static Stencil stencil(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Returns the flux stencil.
Definition: fluxstencil.hh:57
The flux stencil specialized for different discretization schemes.
Definition: fluxstencil.hh:33
The available discretization methods in Dumux.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27