version 3.10-dev
fluxstencil.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_DISCRETIZATION_FLUXSTENCIL_HH
13#define DUMUX_DISCRETIZATION_FLUXSTENCIL_HH
14
15#include <vector>
16
17#include <dune/common/reservedvector.hh>
20
21namespace Dumux {
22
32template<class FVElementGeometry, class DiscretizationMethod = typename FVElementGeometry::GridGeometry::DiscretizationMethod>
34
35/*
36 * \ingroup Discretization
37 * \brief Flux stencil specialization for the cell-centered tpfa scheme
38 * \tparam FVElementGeometry The local view on the finite volume grid geometry
39 */
40template<class FVElementGeometry>
41class FluxStencil<FVElementGeometry, DiscretizationMethods::CCTpfa>
42{
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;
47 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
48
49public:
51 using ScvfStencilIForJ = Dune::ReservedVector<GridIndexType, 1>;
52
54 using Stencil = typename SubControlVolumeFace::Traits::GridIndexStorage;
55
57 static Stencil stencil(const Element& element,
58 const FVElementGeometry& fvGeometry,
59 const SubControlVolumeFace& scvf)
60 {
61 if (scvf.boundary())
62 return Stencil({scvf.insideScvIdx()});
63 else if (scvf.numOutsideScvs() > 1)
64 {
65 Stencil stencil({scvf.insideScvIdx()});
66 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
67 stencil.push_back(scvf.outsideScvIdx(i));
68 return stencil;
69 }
70 else
71 return Stencil({scvf.insideScvIdx(), scvf.outsideScvIdx()});
72 }
73};
74
75/*
76 * \ingroup Discretization
77 * \brief Flux stencil specialization for the cell-centered mpfa scheme
78 * \tparam FVElementGeometry The local view on the finite volume grid geometry
79 */
80template<class FVElementGeometry>
81class FluxStencil<FVElementGeometry, DiscretizationMethods::CCMpfa>
82{
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;
87 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
88
89 // Use the stencil type of the primary interaction volume
90 using NodalIndexSet = typename GridGeometry::GridIVIndexSets::DualGridIndexSet::NodalIndexSet;
91
92public:
94 using ScvfStencilIForJ = std::vector<GridIndexType>;
95
97 using Stencil = typename NodalIndexSet::NodalGridStencilType;
98
100 static const Stencil& stencil(const Element& element,
101 const FVElementGeometry& fvGeometry,
102 const SubControlVolumeFace& scvf)
103 {
104 const auto& gridGeometry = fvGeometry.gridGeometry();
105
106 // return the scv (element) indices in the interaction region
107 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
108 return gridGeometry.gridInteractionVolumeIndexSets().secondaryIndexSet(scvf).gridScvIndices();
109 else
110 return gridGeometry.gridInteractionVolumeIndexSets().primaryIndexSet(scvf).gridScvIndices();
111 }
112};
113
114} // end namespace Dumux
115
116#endif
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
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27