version 3.8
multidomain/facet/cellcentered/mpfa/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_MULTIDOMAIN_FACET_CC_MPFA_O_INTERACTIONVOLUME_HH
13#define DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_INTERACTIONVOLUME_HH
14
20
21#include "localassembler.hh"
23
24namespace Dumux {
25
27template< class Traits > class CCMpfaOFacetCouplingInteractionVolume;
28
38template< class NodalIndexSet, class Scalar >
40: public CCMpfaODefaultInteractionVolumeTraits< NodalIndexSet, Scalar >
41{
42private:
45
46 static constexpr int dim = NodalIndexSet::Traits::GridView::dimension;
47 static constexpr int dimWorld = NodalIndexSet::Traits::GridView::dimensionworld;
48
49public:
54
56 template<class Problem, class FVElementGeometry, class ElemVolVars>
58};
59
66template< class Traits >
68: public CCMpfaInteractionVolumeBase< Traits >
69{
70 using GridView = typename Traits::GridView;
71 using Element = typename GridView::template Codim<0>::Entity;
72
73 using IndexSet = typename Traits::IndexSet;
74 using GridIndexType = typename IndexSet::GridIndexType;
75 using LocalIndexType = typename IndexSet::LocalIndexType;
76 using Stencil = typename IndexSet::NodalGridStencilType;
77
78 using LocalScvType = typename Traits::LocalScvType;
79 using LocalScvfType = typename Traits::LocalScvfType;
80 using LocalFaceData = typename Traits::LocalFaceData;
81
82public:
85
88 {
89 GridIndexType scvfIdx_;
90
91 public:
93 InteriorBoundaryData(GridIndexType scvfIdx) : scvfIdx_(scvfIdx) {}
94
96 GridIndexType scvfIndex() const { return scvfIdx_; }
97 };
98
101
103 template< class Problem, class FVElementGeometry >
104 void bind(const IndexSet& indexSet,
105 const Problem& problem,
106 const FVElementGeometry& fvGeometry)
107 {
108 // for the o-scheme, the stencil is equal to the scv
109 // index set of the dual grid's nodal index set
110 stencil_ = &indexSet.nodalIndexSet().gridScvIndices();
111
112 // find out how many facet elements appear in this iv
113 std::size_t numFacetElems = 0;
114 std::size_t numOutsideFaces = 0;
115 std::vector<bool> isOnInteriorBoundary(indexSet.numFaces(), false);
116 for (LocalIndexType fIdx = 0; fIdx < indexSet.numFaces(); ++fIdx)
117 {
118 const auto& scvf = fvGeometry.scvf(indexSet.gridScvfIndex(fIdx));
119 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
120 if (problem.couplingManager().isOnInteriorBoundary(element, scvf))
121 {
122 numFacetElems++;
123 numOutsideFaces += scvf.numOutsideScvs();
124 isOnInteriorBoundary[fIdx] = true;
125 interiorBoundaryData_.emplace_back( scvf.index() );
126 }
127 }
128
129 // number of interaction-volume-local scvs(=node-local for o-scheme) and scvfs
130 numFaces_ = indexSet.numFaces() + numOutsideFaces;
131 const auto numLocalScvs = indexSet.numScvs();
132 const auto numGlobalScvfs = indexSet.nodalIndexSet().numScvfs();
133
134 // reserve memory for local entities
135 elements_.clear(); elements_.reserve(numLocalScvs);
136 scvs_.clear(); scvs_.reserve(numLocalScvs);
137 scvfs_.clear(); scvfs_.reserve(numFaces_);
138 localFaceData_.clear(); localFaceData_.reserve(numGlobalScvfs);
139 dirichletData_.clear(); dirichletData_.reserve(numFaces_);
140
141 // keep track of the number of unknowns etc
142 numUnknowns_ = 0;
143 numKnowns_ = numLocalScvs + numFacetElems;
144
145 // index map from grid scvf index to local scvf index
146 std::unordered_map<GridIndexType, LocalIndexType> scvfIndexMap;
147
148 // set up objects related to sub-control volume faces
149 LocalIndexType facetElementCounter = 0;
150 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < indexSet.numFaces(); ++faceIdxLocal)
151 {
152 const auto gridScvfIdx = indexSet.gridScvfIndex(faceIdxLocal);
153 const auto& flipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[gridScvfIdx];
154 const auto& scvf = fvGeometry.scvf(gridScvfIdx);
155 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
156
157 // the neighboring scvs in local indices (order: 0 - inside scv, 1..n - outside scvs)
158 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
159 const auto numNeighborScvs = neighborScvIndicesLocal.size();
160
161 // the ív-local scvf index of the face about to be created
162 const auto curLocalScvfIdx = scvfs_.size();
163 scvfIndexMap[gridScvfIdx] = curLocalScvfIdx;
164 localFaceData_.emplace_back(curLocalScvfIdx, neighborScvIndicesLocal[0], scvf.index());
165
166 // on interior boundaries, create local scvfs for inside AND all outside scvfs
167 if (isOnInteriorBoundary[faceIdxLocal])
168 {
169 const LocalIndexType facetLocalDofIdx = numLocalScvs + facetElementCounter++;
170 const bool isDirichlet = problem.interiorBoundaryTypes(element, scvf).hasOnlyDirichlet();
171
172 // create local scvf
173 if (isDirichlet)
174 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, facetLocalDofIdx, /*isDirichlet*/true, facetLocalDofIdx);
175 else
176 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false, facetLocalDofIdx);
177
178 // create "outside" local scvfs
179 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
180 {
181 const auto outsideGridScvfIdx = flipScvfIdxSet[i-1];
182 const auto& flipScvf = fvGeometry.scvf(outsideGridScvfIdx);
183 const auto& outsideFlipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[outsideGridScvfIdx];
184
185 // rearrange the neighbor scv index vector corresponding to this scvfs flip scvf index set
186 using std::swap;
187 auto outsideNeighborScvIdxSet = neighborScvIndicesLocal;
188 outsideNeighborScvIdxSet[0] = outsideNeighborScvIdxSet[i];
189 for (LocalIndexType j = 0; j < outsideFlipScvfIdxSet.size(); ++j)
190 {
191 const auto flipScvfIdx = outsideFlipScvfIdxSet[j];
192 auto it = std::find(flipScvfIdxSet.begin(), flipScvfIdxSet.end(), flipScvfIdx);
193
194 // if we found the index, use corresponding local scv index
195 if (it != flipScvfIdxSet.end())
196 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[std::distance(flipScvfIdxSet.begin(), it)+1];
197
198 // otherwise this must be the "inside" scvf again
199 else
200 {
201 assert(flipScvfIdx == gridScvfIdx);
202 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[0];
203 }
204 }
205
206 scvfIndexMap[outsideGridScvfIdx] = curLocalScvfIdx+i;
207 localFaceData_.emplace_back(curLocalScvfIdx+i, outsideNeighborScvIdxSet[0], flipScvf.index());
208 if (isDirichlet)
209 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, facetLocalDofIdx, /*isDirichlet*/true, facetLocalDofIdx);
210 else
211 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, numUnknowns_++, /*isDirichlet*/false, facetLocalDofIdx);
212 }
213 }
214
215 // otherwise create boundary scvf ...
216 else if (scvf.boundary())
217 {
218 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
219 {
220 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numKnowns_++, /*isDirichlet*/true);
221 dirichletData_.emplace_back(scvf.outsideScvIdx());
222 }
223 else
224 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false);
225 }
226
227 // ... or interior scvf
228 else
229 {
230 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false);
231
232 // add local face data objects for the outside faces
233 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
234 {
235 // loop over scvfs in outside scv until we find the one coinciding with current scvf
236 const auto outsideLocalScvIdx = neighborScvIndicesLocal[i];
237 const auto& flipScvfIndex = fvGeometry.gridGeometry().flipScvfIndexSet()[scvf.index()][i-1];
238 const auto& flipScvf = fvGeometry.scvf(flipScvfIndex);
239 scvfIndexMap[flipScvfIndex] = curLocalScvfIdx;
240 localFaceData_.emplace_back(curLocalScvfIdx, // iv-local scvf idx
241 outsideLocalScvIdx, // iv-local scv index
242 i-1, // scvf-local index in outside faces
243 flipScvf.index()); // global scvf index
244 }
245 }
246 }
247
248 // set up stuff related to sub-control volumes
249 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numLocalScvs; scvIdxLocal++)
250 {
251 elements_.emplace_back(fvGeometry.gridGeometry().element( stencil()[scvIdxLocal] ));
252 scvs_.emplace_back(fvGeometry.gridGeometry().mpfaHelper(),
253 fvGeometry,
254 fvGeometry.scv( stencil()[scvIdxLocal] ),
255 scvIdxLocal,
256 indexSet,
257 scvfIndexMap);
258 }
259 }
260
262 std::size_t numFaces() const
263 { return numFaces_; }
264
266 std::size_t numUnknowns() const
267 { return numUnknowns_; }
268
270 std::size_t numKnowns() const
271 { return numKnowns_; }
272
274 std::size_t numScvs() const
275 { return scvs_.size(); }
276
278 const Stencil& stencil() const
279 { return *stencil_; }
280
282 const Element& element(LocalIndexType ivLocalScvIdx) const
283 { return elements_[ivLocalScvIdx]; }
284
286 const LocalScvfType& localScvf(LocalIndexType ivLocalScvfIdx) const
287 { return scvfs_[ivLocalScvfIdx]; }
288
290 const LocalScvType& localScv(LocalIndexType ivLocalScvIdx) const
291 { return scvs_[ivLocalScvIdx]; }
292
294 const std::vector<LocalFaceData>& localFaceData() const
295 { return localFaceData_; }
296
298 const std::vector<DirichletData>& dirichletData() const
299 { return dirichletData_; }
300
302 const std::vector<InteriorBoundaryData>& interiorBoundaryData() const
303 { return interiorBoundaryData_; }
304
306 template< class FVElementGeometry >
307 auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry& fvGeometry) const
308 { return CCMpfaOScvGeometryHelper<LocalScvType>::computeScvGeometry(ivLocalScvIdx, *this, fvGeometry); }
309
311 template< class NI >
312 static constexpr std::size_t numIVAtVertex(const NI& nodalIndexSet)
313 { return 1; }
314
317 template< class IvIndexSetContainer,
318 class ScvfIndexMap,
319 class NodalIndexSet,
320 class FlipScvfIndexSet >
321 static void addIVIndexSets(IvIndexSetContainer& ivIndexSetContainer,
322 ScvfIndexMap& scvfIndexMap,
323 const NodalIndexSet& nodalIndexSet,
324 const FlipScvfIndexSet& flipScvfIndexSet)
325 {
326 // reuse the function of the standard mpfa-o interaction volume
328 scvfIndexMap,
329 nodalIndexSet,
330 flipScvfIndexSet);
331 }
332
333private:
334 // pointer to cell stencil (in iv index set)
335 const Stencil* stencil_;
336
337 // Variables defining the local scope
338 std::vector<Element> elements_;
339 std::vector<LocalScvType> scvs_;
340 std::vector<LocalScvfType> scvfs_;
341 std::vector<LocalFaceData> localFaceData_;
342 std::vector<DirichletData> dirichletData_;
343 std::vector<InteriorBoundaryData> interiorBoundaryData_;
344
345 // sizes involved in the local system equations
346 std::size_t numFaces_;
347 std::size_t numUnknowns_;
348 std::size_t numKnowns_;
349};
350
351} // end namespace
352
353#endif
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:56
Define data structure to store which scvfs lie on interior boundaries.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:88
GridIndexType scvfIndex() const
Return corresponding scvf index.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:96
InteriorBoundaryData(GridIndexType scvfIdx)
Constructor.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:93
Forward declaration of the facet coupling MPFA-O interaction volume.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:69
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:274
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:286
const std::vector< LocalFaceData > & localFaceData() const
returns a reference to the container with the local face data
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:294
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:307
const std::vector< DirichletData > & dirichletData() const
returns a reference to the information container on Dirichlet BCs within this iv
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:298
typename CCMpfaOInteractionVolume< Traits >::DirichletData DirichletData
Reuse standard o-scheme's Dirichlet Data class.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:84
const std::vector< InteriorBoundaryData > & interiorBoundaryData() const
returns a reference to the data container on interior boundaries
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:302
std::size_t numKnowns() const
returns the number of (in this context) known solution values within this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:270
void bind(const IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Sets up the local scope for a given iv index set.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:104
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the grid element corresponding to a given iv-local scv idx
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:282
static constexpr std::size_t numIVAtVertex(const NI &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:312
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:266
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:278
std::size_t numFaces() const
returns the number of primary scvfs of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:262
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:290
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:321
static constexpr MpfaMethods MpfaMethod
publicly state the mpfa-scheme this interaction volume is associated with
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:100
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme in the context of...
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:35
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:120
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:267
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
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: multidomain/facet/cellcentered/mpfa/localassembler.hh:43
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
Class for the interaction volume of the mpfa-o scheme.
MpfaMethods
The available mpfa schemes in Dumux.
Definition: methods.hh:22
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
Base class for interaction volumes of mpfa methods.
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
The default interaction volume traits class for the mpfa-o method in the context of facet coupling....
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:41
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme in the conte...
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:94
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:126