version 3.10-dev
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 if (!scvf.boundary())
124 numOutsideFaces += scvf.numOutsideScvs();
125 isOnInteriorBoundary[fIdx] = true;
126 interiorBoundaryData_.emplace_back( scvf.index() );
127 }
128 }
129
130 // number of interaction-volume-local scvs(=node-local for o-scheme) and scvfs
131 numFaces_ = indexSet.numFaces() + numOutsideFaces;
132 const auto numLocalScvs = indexSet.numScvs();
133 const auto numGlobalScvfs = indexSet.nodalIndexSet().numScvfs();
134
135 // reserve memory for local entities
136 elements_.clear(); elements_.reserve(numLocalScvs);
137 scvs_.clear(); scvs_.reserve(numLocalScvs);
138 scvfs_.clear(); scvfs_.reserve(numFaces_);
139 localFaceData_.clear(); localFaceData_.reserve(numGlobalScvfs);
140 dirichletData_.clear(); dirichletData_.reserve(numFaces_);
141
142 // keep track of the number of unknowns etc
143 numUnknowns_ = 0;
144 numKnowns_ = numLocalScvs + numFacetElems;
145
146 // index map from grid scvf index to local scvf index
147 std::unordered_map<GridIndexType, LocalIndexType> scvfIndexMap;
148
149 // set up objects related to sub-control volume faces
150 LocalIndexType facetElementCounter = 0;
151 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < indexSet.numFaces(); ++faceIdxLocal)
152 {
153 const auto gridScvfIdx = indexSet.gridScvfIndex(faceIdxLocal);
154 const auto& flipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[gridScvfIdx];
155 const auto& scvf = fvGeometry.scvf(gridScvfIdx);
156 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
157
158 // the neighboring scvs in local indices (order: 0 - inside scv, 1..n - outside scvs)
159 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
160 const auto numNeighborScvs = neighborScvIndicesLocal.size();
161
162 // the ív-local scvf index of the face about to be created
163 const auto curLocalScvfIdx = scvfs_.size();
164 scvfIndexMap[gridScvfIdx] = curLocalScvfIdx;
165 localFaceData_.emplace_back(curLocalScvfIdx, neighborScvIndicesLocal[0], scvf.index());
166
167 // on interior boundaries, create local scvfs for inside AND all outside scvfs
168 if (isOnInteriorBoundary[faceIdxLocal])
169 {
170 const LocalIndexType facetLocalDofIdx = numLocalScvs + facetElementCounter++;
171 const bool isDirichlet = problem.interiorBoundaryTypes(element, scvf).hasOnlyDirichlet();
172
173 // create local scvf
174 if (isDirichlet)
175 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, facetLocalDofIdx, /*isDirichlet*/true, facetLocalDofIdx);
176 else
177 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false, facetLocalDofIdx);
178
179 // create "outside" local scvfs
180 if (!scvf.boundary())
181 {
182 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
183 {
184 const auto outsideGridScvfIdx = flipScvfIdxSet[i-1];
185 const auto& flipScvf = fvGeometry.scvf(outsideGridScvfIdx);
186 const auto& outsideFlipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[outsideGridScvfIdx];
187
188 // rearrange the neighbor scv index vector corresponding to this scvfs flip scvf index set
189 using std::swap;
190 auto outsideNeighborScvIdxSet = neighborScvIndicesLocal;
191 outsideNeighborScvIdxSet[0] = outsideNeighborScvIdxSet[i];
192 for (LocalIndexType j = 0; j < outsideFlipScvfIdxSet.size(); ++j)
193 {
194 const auto flipScvfIdx = outsideFlipScvfIdxSet[j];
195 auto it = std::find(flipScvfIdxSet.begin(), flipScvfIdxSet.end(), flipScvfIdx);
196
197 // if we found the index, use corresponding local scv index
198 if (it != flipScvfIdxSet.end())
199 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[std::distance(flipScvfIdxSet.begin(), it)+1];
200
201 // otherwise this must be the "inside" scvf again
202 else
203 {
204 assert(flipScvfIdx == gridScvfIdx);
205 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[0];
206 }
207 }
208
209 scvfIndexMap[outsideGridScvfIdx] = curLocalScvfIdx+i;
210 localFaceData_.emplace_back(curLocalScvfIdx+i, outsideNeighborScvIdxSet[0], flipScvf.index());
211 if (isDirichlet)
212 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, facetLocalDofIdx, /*isDirichlet*/true, facetLocalDofIdx);
213 else
214 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, numUnknowns_++, /*isDirichlet*/false, facetLocalDofIdx);
215 }
216 }
217 }
218
219 // otherwise create boundary scvf ...
220 else if (scvf.boundary())
221 {
222 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
223 {
224 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numKnowns_++, /*isDirichlet*/true);
225 dirichletData_.emplace_back(scvf.outsideScvIdx());
226 }
227 else
228 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false);
229 }
230
231 // ... or interior scvf
232 else
233 {
234 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false);
235
236 // add local face data objects for the outside faces
237 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
238 {
239 // loop over scvfs in outside scv until we find the one coinciding with current scvf
240 const auto outsideLocalScvIdx = neighborScvIndicesLocal[i];
241 const auto& flipScvfIndex = fvGeometry.gridGeometry().flipScvfIndexSet()[scvf.index()][i-1];
242 const auto& flipScvf = fvGeometry.scvf(flipScvfIndex);
243 scvfIndexMap[flipScvfIndex] = curLocalScvfIdx;
244 localFaceData_.emplace_back(curLocalScvfIdx, // iv-local scvf idx
245 outsideLocalScvIdx, // iv-local scv index
246 i-1, // scvf-local index in outside faces
247 flipScvf.index()); // global scvf index
248 }
249 }
250 }
251
252 // set up stuff related to sub-control volumes
253 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numLocalScvs; scvIdxLocal++)
254 {
255 elements_.emplace_back(fvGeometry.gridGeometry().element( stencil()[scvIdxLocal] ));
256 scvs_.emplace_back(fvGeometry.gridGeometry().mpfaHelper(),
257 fvGeometry,
258 fvGeometry.scv( stencil()[scvIdxLocal] ),
259 scvIdxLocal,
260 indexSet,
261 scvfIndexMap);
262 }
263 }
264
266 std::size_t numFaces() const
267 { return numFaces_; }
268
270 std::size_t numUnknowns() const
271 { return numUnknowns_; }
272
274 std::size_t numKnowns() const
275 { return numKnowns_; }
276
278 std::size_t numScvs() const
279 { return scvs_.size(); }
280
282 const Stencil& stencil() const
283 { return *stencil_; }
284
286 const Element& element(LocalIndexType ivLocalScvIdx) const
287 { return elements_[ivLocalScvIdx]; }
288
290 const LocalScvfType& localScvf(LocalIndexType ivLocalScvfIdx) const
291 { return scvfs_[ivLocalScvfIdx]; }
292
294 const LocalScvType& localScv(LocalIndexType ivLocalScvIdx) const
295 { return scvs_[ivLocalScvIdx]; }
296
298 const std::vector<LocalFaceData>& localFaceData() const
299 { return localFaceData_; }
300
302 const std::vector<DirichletData>& dirichletData() const
303 { return dirichletData_; }
304
306 const std::vector<InteriorBoundaryData>& interiorBoundaryData() const
307 { return interiorBoundaryData_; }
308
310 template< class FVElementGeometry >
311 auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry& fvGeometry) const
312 { return CCMpfaOScvGeometryHelper<LocalScvType>::computeScvGeometry(ivLocalScvIdx, *this, fvGeometry); }
313
315 template< class NI >
316 static constexpr std::size_t numIVAtVertex(const NI& nodalIndexSet)
317 { return 1; }
318
321 template< class IvIndexSetContainer,
322 class ScvfIndexMap,
323 class NodalIndexSet,
324 class FlipScvfIndexSet >
325 static void addIVIndexSets(IvIndexSetContainer& ivIndexSetContainer,
326 ScvfIndexMap& scvfIndexMap,
327 const NodalIndexSet& nodalIndexSet,
328 const FlipScvfIndexSet& flipScvfIndexSet)
329 {
330 // reuse the function of the standard mpfa-o interaction volume
332 scvfIndexMap,
333 nodalIndexSet,
334 flipScvfIndexSet);
335 }
336
337private:
338 // pointer to cell stencil (in iv index set)
339 const Stencil* stencil_;
340
341 // Variables defining the local scope
342 std::vector<Element> elements_;
343 std::vector<LocalScvType> scvs_;
344 std::vector<LocalScvfType> scvfs_;
345 std::vector<LocalFaceData> localFaceData_;
346 std::vector<DirichletData> dirichletData_;
347 std::vector<InteriorBoundaryData> interiorBoundaryData_;
348
349 // sizes involved in the local system equations
350 std::size_t numFaces_;
351 std::size_t numUnknowns_;
352 std::size_t numKnowns_;
353};
354
355} // end namespace
356
357#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:278
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:290
const std::vector< LocalFaceData > & localFaceData() const
returns a reference to the container with the local face data
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:298
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:311
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:302
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:306
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:274
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:286
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:316
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:270
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:282
std::size_t numFaces() const
returns the number of primary scvfs of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:266
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:294
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:325
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