version 3.10-dev
discretization/cellcentered/mpfa/omethod/interactionvolume.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_CC_MPFA_O_INTERACTIONVOLUME_HH
13#define DUMUX_DISCRETIZATION_CC_MPFA_O_INTERACTIONVOLUME_HH
14
15#include <type_traits>
16
17#include <dune/common/dynmatrix.hh>
18#include <dune/common/dynvector.hh>
19#include <dune/common/fvector.hh>
20#include <dune/common/reservedvector.hh>
21
25
26#include "localassembler.hh"
29#include "scvgeometryhelper.hh"
30
31namespace Dumux {
32
34template< class Traits > class CCMpfaOInteractionVolume;
35
46template< class NodalIndexSet, class Scalar >
48{
49private:
50 using GridIndexType = typename NodalIndexSet::GridIndexType;
51 using LocalIndexType = typename NodalIndexSet::LocalIndexType;
52
53 static constexpr int dim = NodalIndexSet::Traits::GridView::dimension;
54 static constexpr int dimWorld = NodalIndexSet::Traits::GridView::dimensionworld;
55
56 using DimVector = Dune::FieldVector<Scalar, dim>;
57 using FaceOmegas = typename std::conditional< (dim<dimWorld),
58 std::vector<DimVector>,
59 Dune::ReservedVector<DimVector, 2> >::type;
60
62 struct MVTraits
63 {
64 using OmegaStorage = std::vector< FaceOmegas >;
65
66 using AMatrix = Dune::DynamicMatrix< Scalar >;
67 using BMatrix = Dune::DynamicMatrix< Scalar >;
68 using CMatrix = Dune::DynamicMatrix< Scalar >;
69 using DMatrix = Dune::DynamicMatrix< Scalar >;
70 using TMatrix = Dune::DynamicMatrix< Scalar >;
71 using CellVector = Dune::DynamicVector< Scalar >;
72 using FaceVector = Dune::DynamicVector< Scalar >;
73 };
74
75public:
77 using GridView = typename NodalIndexSet::Traits::GridView;
87 using MatVecTraits = MVTraits;
88
90 template<class Problem, class FVElementGeometry, class ElemVolVars>
92};
93
100template< class Traits >
102: public CCMpfaInteractionVolumeBase< Traits >
103{
104 using GridView = typename Traits::GridView;
105 using Element = typename GridView::template Codim<0>::Entity;
106
107 using IndexSet = typename Traits::IndexSet;
108 using GridIndexType = typename IndexSet::GridIndexType;
109 using LocalIndexType = typename IndexSet::LocalIndexType;
110 using Stencil = typename IndexSet::NodalGridStencilType;
111
112 using LocalScvType = typename Traits::LocalScvType;
113 using LocalScvfType = typename Traits::LocalScvfType;
114 using LocalFaceData = typename Traits::LocalFaceData;
115
116public:
120 {
121 GridIndexType volVarIndex_;
122 public:
124 DirichletData(const GridIndexType index) : volVarIndex_(index) {}
125
127 GridIndexType volVarIndex() const { return volVarIndex_; }
128 };
129
132
134 template< class Problem, class FVElementGeometry >
135 void bind(const IndexSet& indexSet,
136 const Problem& problem,
137 const FVElementGeometry& fvGeometry)
138 {
139 // for the o-scheme, the stencil is equal to the scv
140 // index set of the dual grid's nodal index set
141 stencil_ = &indexSet.nodalIndexSet().gridScvIndices();
142
143 // number of interaction-volume-local scvs(=node-local for o-scheme) and scvfs
144 numFaces_ = indexSet.numFaces();
145 const auto numLocalScvs = indexSet.numScvs();
146 const auto numGlobalScvfs = indexSet.nodalIndexSet().numScvfs();
147
148 // reserve memory for local entities
149 elements_.clear(); elements_.reserve(numLocalScvs);
150 scvs_.clear(); scvs_.reserve(numLocalScvs);
151 scvfs_.clear(); scvfs_.reserve(numFaces_);
152 localFaceData_.clear(); localFaceData_.reserve(numGlobalScvfs);
153 dirichletData_.clear(); dirichletData_.reserve(numFaces_);
154
155 // set up stuff related to sub-control volumes
156 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numLocalScvs; scvIdxLocal++)
157 {
158 elements_.emplace_back(fvGeometry.gridGeometry().element( stencil()[scvIdxLocal] ));
159 scvs_.emplace_back(fvGeometry.gridGeometry().mpfaHelper(),
160 fvGeometry,
161 fvGeometry.scv( stencil()[scvIdxLocal] ),
162 scvIdxLocal,
163 indexSet);
164 }
165
166 // keep track of the number of unknowns etc
167 numUnknowns_ = 0;
168 numKnowns_ = numLocalScvs;
169
170 // set up quantities related to sub-control volume faces
171 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < numFaces_; ++faceIdxLocal)
172 {
173 const auto& scvf = fvGeometry.scvf(indexSet.gridScvfIndex(faceIdxLocal));
174
175 // the neighboring scvs in local indices (order: 0 - inside scv, 1..n - outside scvs)
176 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
177 const auto numNeighborScvs = neighborScvIndicesLocal.size();
178 localFaceData_.emplace_back(faceIdxLocal, neighborScvIndicesLocal[0], scvf.index());
179
180 // create iv-local scvf object
181 if (scvf.boundary())
182 {
183 if (problem.boundaryTypes(elements_[neighborScvIndicesLocal[0]], scvf).hasOnlyDirichlet())
184 {
185 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numKnowns_++, /*isDirichlet*/true);
186 dirichletData_.emplace_back(scvf.outsideScvIdx());
187 }
188 else
189 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false);
190 }
191 else
192 {
193 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false);
194
195 // add local face data objects for the outside faces
196 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
197 {
198 // loop over scvfs in outside scv until we find the one coinciding with current scvf
199 const auto outsideLocalScvIdx = neighborScvIndicesLocal[i];
200 const auto& flipScvfIndex = fvGeometry.gridGeometry().flipScvfIndexSet()[scvf.index()][i-1];
201 const auto& flipScvf = fvGeometry.scvf(flipScvfIndex);
202 localFaceData_.emplace_back(faceIdxLocal, // iv-local scvf idx
203 outsideLocalScvIdx, // iv-local scv index
204 i-1, // scvf-local index in outside faces
205 flipScvf.index()); // global scvf index
206 }
207 }
208 }
209 }
210
212 std::size_t numFaces() const
213 { return numFaces_; }
214
216 std::size_t numUnknowns() const
217 { return numUnknowns_; }
218
220 std::size_t numKnowns() const
221 { return numKnowns_; }
222
224 std::size_t numScvs() const
225 { return scvs_.size(); }
226
228 const Stencil& stencil() const
229 { return *stencil_; }
230
232 const Element& element(LocalIndexType ivLocalScvIdx) const
233 { return elements_[ivLocalScvIdx]; }
234
236 const LocalScvfType& localScvf(LocalIndexType ivLocalScvfIdx) const
237 { return scvfs_[ivLocalScvfIdx]; }
238
240 const LocalScvType& localScv(LocalIndexType ivLocalScvIdx) const
241 { return scvs_[ivLocalScvIdx]; }
242
244 const std::vector<LocalFaceData>& localFaceData() const
245 { return localFaceData_; }
246
248 const std::vector<DirichletData>& dirichletData() const
249 { return dirichletData_; }
250
252 template< class FVElementGeometry >
253 auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry& fvGeometry) const
254 { return CCMpfaOScvGeometryHelper<LocalScvType>::computeScvGeometry(ivLocalScvIdx, *this, fvGeometry); }
255
257 template< class NI >
258 static constexpr std::size_t numIVAtVertex(const NI& nodalIndexSet)
259 { return 1; }
260
263 template< class IvIndexSetContainer,
264 class ScvfIndexMap,
265 class NodalIndexSet,
266 class FlipScvfIndexSet >
267 static void addIVIndexSets(IvIndexSetContainer& ivIndexSetContainer,
268 ScvfIndexMap& scvfIndexMap,
269 const NodalIndexSet& nodalIndexSet,
270 const FlipScvfIndexSet& flipScvfIndexSet)
271 {
272 // the global index of the iv index set that is about to be created
273 const auto curGlobalIndex = ivIndexSetContainer.size();
274
275 // make the one index set for this node
276 ivIndexSetContainer.emplace_back(nodalIndexSet, flipScvfIndexSet);
277
278 // store the index mapping
279 for (const auto scvfIdx : nodalIndexSet.gridScvfIndices())
280 scvfIndexMap[scvfIdx] = curGlobalIndex;
281 }
282
283private:
284 // pointer to cell stencil (in iv index set)
285 const Stencil* stencil_;
286
287 // Variables defining the local scope
288 std::vector<Element> elements_;
289 std::vector<LocalScvType> scvs_;
290 std::vector<LocalScvfType> scvfs_;
291 std::vector<LocalFaceData> localFaceData_;
292 std::vector<DirichletData> dirichletData_;
293
294 // sizes involved in the local system equations
295 std::size_t numFaces_;
296 std::size_t numUnknowns_;
297 std::size_t numKnowns_;
298};
299
300} // end namespace Dumux
301
302#endif
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:56
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:120
GridIndexType volVarIndex() const
Return corresponding vol var index.
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:127
DirichletData(const GridIndexType index)
Constructor.
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:124
Forward declaration of the o-method's interaction volume.
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:103
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:236
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:216
static constexpr std::size_t numIVAtVertex(const NI &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:258
void bind(const IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Sets up the local scope for a given iv index set.
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:135
static constexpr MpfaMethods MpfaMethod
publicly state the mpfa-scheme this interaction volume is associated with
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:131
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:267
std::size_t numKnowns() const
returns the number of (in this context) known solution values within this interaction volume
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:220
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the grid element corresponding to a given iv-local scv idx
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:232
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:240
const std::vector< LocalFaceData > & localFaceData() const
returns a reference to the container with the local face data
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:244
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:228
const std::vector< DirichletData > & dirichletData() const
returns a reference to the information container on Dirichlet BCs within this iv
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:248
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:224
std::size_t numFaces() const
returns the number of primary scvfs of this interaction volume
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:212
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:253
The interaction volume index set class for the mpfa-o scheme.
Definition: interactionvolumeindexset.hh:29
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:31
static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const InteractionVolume &iv, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: scvgeometryhelper.hh:61
General implementation of a data structure holding interaction volume-local information for a grid su...
Definition: localfacedata.hh:30
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: discretization/cellcentered/mpfa/omethod/localassembler.hh:36
MpfaMethods
The available mpfa schemes in Dumux.
Definition: methods.hh:22
Base class for interaction volumes of mpfa methods.
Class for the index set within an interaction volume of the mpfa-o scheme.
Data structure holding interaction volume-local information for a grid subb-control volume face embed...
The available mpfa schemes in Dumux.
Classes for sub control entities of the mpfa-o method in the context of facet coupling.
Definition: adapt.hh:17
Helper class providing functionality to compute the geometry of the interaction-volume local sub-cont...
The default interaction volume traits class for the mpfa-o method. This uses dynamic types types for ...
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:48
MVTraits MatVecTraits
export the matrix/vector traits to be used by the iv
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:87
typename NodalIndexSet::Traits::GridView GridView
export the type of grid view
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:77
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:126