3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 2 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_INTERACTIONVOLUME_HH
25#define DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_INTERACTIONVOLUME_HH
26
32
33#include "localassembler.hh"
35
36namespace Dumux {
37
39template< class Traits > class CCMpfaOFacetCouplingInteractionVolume;
40
50template< class NodalIndexSet, class Scalar >
52: public CCMpfaODefaultInteractionVolumeTraits< NodalIndexSet, Scalar >
53{
54private:
57
58 static constexpr int dim = NodalIndexSet::Traits::GridView::dimension;
59 static constexpr int dimWorld = NodalIndexSet::Traits::GridView::dimensionworld;
60
61public:
66
68 template<class Problem, class FVElementGeometry, class ElemVolVars>
70};
71
78template< class Traits >
80: public CCMpfaInteractionVolumeBase< Traits >
81{
82 using GridView = typename Traits::GridView;
83 using Element = typename GridView::template Codim<0>::Entity;
84
85 using IndexSet = typename Traits::IndexSet;
86 using GridIndexType = typename IndexSet::GridIndexType;
87 using LocalIndexType = typename IndexSet::LocalIndexType;
88 using Stencil = typename IndexSet::NodalGridStencilType;
89
90 using LocalScvType = typename Traits::LocalScvType;
91 using LocalScvfType = typename Traits::LocalScvfType;
92 using LocalFaceData = typename Traits::LocalFaceData;
93
94public:
97
100 {
101 GridIndexType scvfIdx_;
102
103 public:
105 InteriorBoundaryData(GridIndexType scvfIdx) : scvfIdx_(scvfIdx) {}
106
108 GridIndexType scvfIndex() const { return scvfIdx_; }
109 };
110
113
115 template< class Problem, class FVElementGeometry >
116 void bind(const IndexSet& indexSet,
117 const Problem& problem,
118 const FVElementGeometry& fvGeometry)
119 {
120 // for the o-scheme, the stencil is equal to the scv
121 // index set of the dual grid's nodal index set
122 stencil_ = &indexSet.nodalIndexSet().gridScvIndices();
123
124 // find out how many facet elements appear in this iv
125 std::size_t numFacetElems = 0;
126 std::size_t numOutsideFaces = 0;
127 std::vector<bool> isOnInteriorBoundary(indexSet.numFaces(), false);
128 for (LocalIndexType fIdx = 0; fIdx < indexSet.numFaces(); ++fIdx)
129 {
130 const auto& scvf = fvGeometry.scvf(indexSet.gridScvfIndex(fIdx));
131 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
132 if (problem.couplingManager().isOnInteriorBoundary(element, scvf))
133 {
134 numFacetElems++;
135 numOutsideFaces += scvf.numOutsideScvs();
136 isOnInteriorBoundary[fIdx] = true;
137 interiorBoundaryData_.emplace_back( scvf.index() );
138 }
139 }
140
141 // number of interaction-volume-local scvs(=node-local for o-scheme) and scvfs
142 numFaces_ = indexSet.numFaces() + numOutsideFaces;
143 const auto numLocalScvs = indexSet.numScvs();
144 const auto numGlobalScvfs = indexSet.nodalIndexSet().numScvfs();
145
146 // reserve memory for local entities
147 elements_.clear(); elements_.reserve(numLocalScvs);
148 scvs_.clear(); scvs_.reserve(numLocalScvs);
149 scvfs_.clear(); scvfs_.reserve(numFaces_);
150 localFaceData_.clear(); localFaceData_.reserve(numGlobalScvfs);
151 dirichletData_.clear(); dirichletData_.reserve(numFaces_);
152
153 // keep track of the number of unknowns etc
154 numUnknowns_ = 0;
155 numKnowns_ = numLocalScvs + numFacetElems;
156
157 // index map from grid scvf index to local scvf index
158 std::unordered_map<GridIndexType, LocalIndexType> scvfIndexMap;
159
160 // set up objects related to sub-control volume faces
161 LocalIndexType facetElementCounter = 0;
162 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < indexSet.numFaces(); ++faceIdxLocal)
163 {
164 const auto gridScvfIdx = indexSet.gridScvfIndex(faceIdxLocal);
165 const auto& flipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[gridScvfIdx];
166 const auto& scvf = fvGeometry.scvf(gridScvfIdx);
167 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
168
169 // the neighboring scvs in local indices (order: 0 - inside scv, 1..n - outside scvs)
170 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
171 const auto numNeighborScvs = neighborScvIndicesLocal.size();
172
173 // the ív-local scvf index of the face about to be created
174 const auto curLocalScvfIdx = scvfs_.size();
175 scvfIndexMap[gridScvfIdx] = curLocalScvfIdx;
176 localFaceData_.emplace_back(curLocalScvfIdx, neighborScvIndicesLocal[0], scvf.index());
177
178 // on interior boundaries, create local scvfs for inside AND all outside scvfs
179 if (isOnInteriorBoundary[faceIdxLocal])
180 {
181 const LocalIndexType facetLocalDofIdx = numLocalScvs + facetElementCounter++;
182 const bool isDirichlet = problem.interiorBoundaryTypes(element, scvf).hasOnlyDirichlet();
183
184 // create local scvf
185 if (isDirichlet)
186 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, facetLocalDofIdx, /*isDirichlet*/true, facetLocalDofIdx);
187 else
188 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false, facetLocalDofIdx);
189
190 // create "outside" local scvfs
191 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
192 {
193 const auto outsideGridScvfIdx = flipScvfIdxSet[i-1];
194 const auto& flipScvf = fvGeometry.scvf(outsideGridScvfIdx);
195 const auto& outsideFlipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[outsideGridScvfIdx];
196
197 // rearrange the neighbor scv index vector corresponding to this scvfs flip scvf index set
198 using std::swap;
199 auto outsideNeighborScvIdxSet = neighborScvIndicesLocal;
200 outsideNeighborScvIdxSet[0] = outsideNeighborScvIdxSet[i];
201 for (LocalIndexType j = 0; j < outsideFlipScvfIdxSet.size(); ++j)
202 {
203 const auto flipScvfIdx = outsideFlipScvfIdxSet[j];
204 auto it = std::find(flipScvfIdxSet.begin(), flipScvfIdxSet.end(), flipScvfIdx);
205
206 // if we found the index, use corresponding local scv index
207 if (it != flipScvfIdxSet.end())
208 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[std::distance(flipScvfIdxSet.begin(), it)+1];
209
210 // otherwise this must be the "inside" scvf again
211 else
212 {
213 assert(flipScvfIdx == gridScvfIdx);
214 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[0];
215 }
216 }
217
218 scvfIndexMap[outsideGridScvfIdx] = curLocalScvfIdx+i;
219 localFaceData_.emplace_back(curLocalScvfIdx+i, outsideNeighborScvIdxSet[0], flipScvf.index());
220 if (isDirichlet)
221 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, facetLocalDofIdx, /*isDirichlet*/true, facetLocalDofIdx);
222 else
223 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, numUnknowns_++, /*isDirichlet*/false, facetLocalDofIdx);
224 }
225 }
226
227 // otherwise crate boundary scvf ...
228 else if (scvf.boundary())
229 {
230 if (problem.boundaryTypes(element, scvf).hasOnlyDirichlet())
231 {
232 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numKnowns_++, /*isDirichlet*/true);
233 dirichletData_.emplace_back(scvf.outsideScvIdx());
234 }
235 else
236 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false);
237 }
238
239 // ... or interior scvf
240 else
241 {
242 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++, /*isDirichlet*/false);
243
244 // add local face data objects for the outside faces
245 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
246 {
247 // loop over scvfs in outside scv until we find the one coinciding with current scvf
248 const auto outsideLocalScvIdx = neighborScvIndicesLocal[i];
249 const auto& flipScvfIndex = fvGeometry.gridGeometry().flipScvfIndexSet()[scvf.index()][i-1];
250 const auto& flipScvf = fvGeometry.scvf(flipScvfIndex);
251 scvfIndexMap[flipScvfIndex] = curLocalScvfIdx;
252 localFaceData_.emplace_back(curLocalScvfIdx, // iv-local scvf idx
253 outsideLocalScvIdx, // iv-local scv index
254 i-1, // scvf-local index in outside faces
255 flipScvf.index()); // global scvf index
256 }
257 }
258 }
259
260 // set up stuff related to sub-control volumes
261 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numLocalScvs; scvIdxLocal++)
262 {
263 elements_.emplace_back(fvGeometry.gridGeometry().element( stencil()[scvIdxLocal] ));
264 scvs_.emplace_back(fvGeometry.gridGeometry().mpfaHelper(),
265 fvGeometry,
266 fvGeometry.scv( stencil()[scvIdxLocal] ),
267 scvIdxLocal,
268 indexSet,
269 scvfIndexMap);
270 }
271 }
272
274 std::size_t numFaces() const
275 { return numFaces_; }
276
278 std::size_t numUnknowns() const
279 { return numUnknowns_; }
280
282 std::size_t numKnowns() const
283 { return numKnowns_; }
284
286 std::size_t numScvs() const
287 { return scvs_.size(); }
288
290 const Stencil& stencil() const
291 { return *stencil_; }
292
294 const Element& element(LocalIndexType ivLocalScvIdx) const
295 { return elements_[ivLocalScvIdx]; }
296
298 const LocalScvfType& localScvf(LocalIndexType ivLocalScvfIdx) const
299 { return scvfs_[ivLocalScvfIdx]; }
300
302 const LocalScvType& localScv(LocalIndexType ivLocalScvIdx) const
303 { return scvs_[ivLocalScvIdx]; }
304
306 const std::vector<LocalFaceData>& localFaceData() const
307 { return localFaceData_; }
308
310 const std::vector<DirichletData>& dirichletData() const
311 { return dirichletData_; }
312
314 const std::vector<InteriorBoundaryData>& interiorBoundaryData() const
315 { return interiorBoundaryData_; }
316
318 template< class FVElementGeometry >
319 auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry& fvGeometry) const
320 { return CCMpfaOScvGeometryHelper<LocalScvType>::computeScvGeometry(ivLocalScvIdx, *this, fvGeometry); }
321
323 template< class NI >
324 static constexpr std::size_t numIVAtVertex(const NI& nodalIndexSet)
325 { return 1; }
326
329 template< class IvIndexSetContainer,
330 class ScvfIndexMap,
331 class NodalIndexSet,
332 class FlipScvfIndexSet >
333 static void addIVIndexSets(IvIndexSetContainer& ivIndexSetContainer,
334 ScvfIndexMap& scvfIndexMap,
335 const NodalIndexSet& nodalIndexSet,
336 const FlipScvfIndexSet& flipScvfIndexSet)
337 {
338 // reuse the function of the standard mpfa-o interaction volume
340 scvfIndexMap,
341 nodalIndexSet,
342 flipScvfIndexSet);
343 }
344
345private:
346 // pointer to cell stencil (in iv index set)
347 const Stencil* stencil_;
348
349 // Variables defining the local scope
350 std::vector<Element> elements_;
351 std::vector<LocalScvType> scvs_;
352 std::vector<LocalScvfType> scvfs_;
353 std::vector<LocalFaceData> localFaceData_;
354 std::vector<DirichletData> dirichletData_;
355 std::vector<InteriorBoundaryData> interiorBoundaryData_;
356
357 // sizes involved in the local system equations
358 std::size_t numFaces_;
359 std::size_t numUnknowns_;
360 std::size_t numKnowns_;
361};
362
363} // end namespace
364
365#endif
Helper class providing functionality to compute the geometry of the interaction-volume local sub-cont...
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.
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:292
MpfaMethods
The available mpfa schemes in Dumux.
Definition: methods.hh:34
Definition: adapt.hh:29
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:68
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:279
The default interaction volume traits class for the mpfa-o method. This uses dynamic types types for ...
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:60
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:132
The interaction volume index set class for the mpfa-o scheme.
Definition: interactionvolumeindexset.hh:41
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: discretization/cellcentered/mpfa/omethod/localassembler.hh:48
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:43
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:138
static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const InteractionVolume &iv, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: scvgeometryhelper.hh:73
Forward declaration of the facet coupling MPFA-O interaction volume.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:81
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:286
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:298
const std::vector< LocalFaceData > & localFaceData() const
returns a reference to the container with the local face data
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:306
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:319
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:310
typename CCMpfaOInteractionVolume< Traits >::DirichletData DirichletData
Reuse standard o-scheme's Dirichlet Data class.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:96
const std::vector< InteriorBoundaryData > & interiorBoundaryData() const
returns a reference to the data container on interior boundaries
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:314
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:282
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:116
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:294
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:324
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:278
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:290
std::size_t numFaces() const
returns the number of primary scvfs of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:274
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:302
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:333
static constexpr MpfaMethods MpfaMethod
publicly state the mpfa-scheme this interaction volume is associated with
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:112
The default interaction volume traits class for the mpfa-o method in the context of facet coupling....
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:53
Define data structure to store which scvfs lie on interior boundaries.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:100
GridIndexType scvfIndex() const
Return corresponding scvf index.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:108
InteriorBoundaryData(GridIndexType scvfIdx)
Constructor.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:105
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: multidomain/facet/cellcentered/mpfa/localassembler.hh:55
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:47
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:106
Class for the interaction volume of the mpfa-o scheme.
Classes for sub control entities of the mpfa-o method in the context of facet coupling.