3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
staticinteractionvolume.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 3 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 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_STATIC_INTERACTIONVOLUME_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_O_STATIC_INTERACTIONVOLUME_HH
28
29#include <dune/common/fmatrix.hh>
30#include <dune/common/fvector.hh>
31#include <dune/common/exceptions.hh>
32
33#include <dumux/common/math.hh>
34
39
40#include "localassembler.hh"
43#include "interactionvolume.hh"
44#include "scvgeometryhelper.hh"
45
46namespace Dumux {
47
49template< class Traits > class CCMpfaOStaticInteractionVolume;
50
63template< class NI, class S, int C, int F >
65{
66private:
67 using GridIndexType = typename NI::GridIndexType;
68 using LocalIndexType = typename NI::LocalIndexType;
69
70 static constexpr int dim = NI::Traits::GridView::dimension;
71 static constexpr int dimWorld = NI::Traits::GridView::dimensionworld;
72
73 using DimVector = Dune::FieldVector<S, dim>;
74 using FaceOmegas = Dune::ReservedVector<DimVector, 2>;
75
77 struct MVTraits
78 {
79 using OmegaStorage = std::array< FaceOmegas, F >;
80
81 using AMatrix = Dune::FieldMatrix< S, F, F >;
82 using BMatrix = Dune::FieldMatrix< S, F, C >;
83 using CMatrix = Dune::FieldMatrix< S, F, F >;
84 using DMatrix = Dune::FieldMatrix< S, F, C >;
85 using TMatrix = Dune::FieldMatrix< S, F, C >;
86 using CellVector = Dune::FieldVector< S, C >;
87 using FaceVector = Dune::FieldVector< S, F >;
88 };
89
90public:
92 using GridView = typename NI::Traits::GridView;
102 using MatVecTraits = MVTraits;
104 static constexpr int numScvs = C;
106 static constexpr int numScvfs = F;
107
109 template<class Problem, class FVElementGeometry, class ElemVolVars>
111};
112
123template< class Traits >
125: public CCMpfaInteractionVolumeBase< Traits >
126{
127 using GridView = typename Traits::GridView;
128 using Element = typename GridView::template Codim<0>::Entity;
129
130 using IndexSet = typename Traits::IndexSet;
131 using GridIndexType = typename IndexSet::GridIndexType;
132 using LocalIndexType = typename IndexSet::LocalIndexType;
133 using Stencil = typename IndexSet::NodalGridStencilType;
134
135 using LocalScvType = typename Traits::LocalScvType;
136 using LocalScvfType = typename Traits::LocalScvfType;
137 using LocalFaceData = typename Traits::LocalFaceData;
138
139 static constexpr int numScvf = Traits::numScvfs;
140 static constexpr int numScv = Traits::numScvs;
141
142public:
144 static_assert(int(GridView::dimension)==int(GridView::dimensionworld), "static iv does not work on surface grids");
145
149 {
150 GridIndexType volVarIndex() const
151 { DUNE_THROW(Dune::InvalidStateException, "Static interaction volume cannot be used on boundaries!"); }
152 };
153
156
158 template< class Problem, class FVElementGeometry >
159 void bind(const IndexSet& indexSet,
160 const Problem& problem,
161 const FVElementGeometry& fvGeometry)
162 {
163 // for the o-scheme, the stencil is equal to the scv
164 // index set of the dual grid's nodal index set
165 assert(indexSet.numScvs() == numScv);
166 stencil_ = &indexSet.nodalIndexSet().gridScvIndices();
167
168 // set up stuff related to sub-control volumes
169 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numScv; scvIdxLocal++)
170 {
171 elements_[scvIdxLocal] = fvGeometry.gridGeometry().element( stencil()[scvIdxLocal] );
172 scvs_[scvIdxLocal] = LocalScvType(fvGeometry.gridGeometry().mpfaHelper(),
173 fvGeometry,
174 fvGeometry.scv( stencil()[scvIdxLocal] ),
175 scvIdxLocal,
176 indexSet);
177 }
178
179 // set up quantities related to sub-control volume faces
180 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < numScvf; ++faceIdxLocal)
181 {
182 const auto& scvf = fvGeometry.scvf(indexSet.gridScvfIndex(faceIdxLocal));
183 assert(!scvf.boundary());
184
185 // the neighboring scvs in local indices (order: 0 - inside scv, 1..n - outside scvs)
186 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
187
188 // create iv-local scvf objects
189 scvfs_[faceIdxLocal] = LocalScvfType(scvf, neighborScvIndicesLocal, faceIdxLocal, /*isDirichlet*/false);
190 localFaceData_[faceIdxLocal*2] = LocalFaceData(faceIdxLocal, neighborScvIndicesLocal[0], scvf.index());
191
192 // add local face data objects for the outside face
193 const auto outsideLocalScvIdx = neighborScvIndicesLocal[1];
194 for (int coord = 0; coord < GridView::dimension; ++coord)
195 {
196 if (indexSet.localScvfIndex(outsideLocalScvIdx, coord) == faceIdxLocal)
197 {
198 const auto globalScvfIdx = indexSet.nodalIndexSet().gridScvfIndex(outsideLocalScvIdx, coord);
199 const auto& flipScvf = fvGeometry.scvf(globalScvfIdx);
200 localFaceData_[faceIdxLocal*2+1] = LocalFaceData(faceIdxLocal, // iv-local scvf idx
201 outsideLocalScvIdx, // iv-local scv index
202 0, // scvf-local index in outside faces
203 flipScvf.index()); // global scvf index
204 break; // go to next outside face
205 }
206 }
207
208 // make sure we found it
209 assert(localFaceData_[faceIdxLocal*2+1].ivLocalInsideScvIndex() == outsideLocalScvIdx);
210 }
211 }
212
214 static constexpr std::size_t numFaces()
215 { return numScvf; }
216
218 static constexpr std::size_t numUnknowns()
219 { return numScvf; }
220
222 static constexpr std::size_t numKnowns()
223 { return numScv; }
224
226 static constexpr std::size_t numScvs()
227 { return numScv; }
228
230 const Stencil& stencil() const
231 { return *stencil_; }
232
234 const Element& element(LocalIndexType ivLocalScvIdx) const
235 { return elements_[ivLocalScvIdx]; }
236
238 const LocalScvfType& localScvf(LocalIndexType ivLocalScvfIdx) const
239 { return scvfs_[ivLocalScvfIdx]; }
240
242 const LocalScvType& localScv(LocalIndexType ivLocalScvIdx) const
243 { return scvs_[ivLocalScvIdx]; }
244
246 const std::array<LocalFaceData, numScvf*2>& localFaceData() const
247 { return localFaceData_; }
248
251 const std::array<DirichletData, 0>& dirichletData() const
252 { return dirichletData_; }
253
255 template< class FVElementGeometry >
256 auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry& fvGeometry) const
257 { return CCMpfaOScvGeometryHelper<LocalScvType>::computeScvGeometry(ivLocalScvIdx, *this, fvGeometry); }
258
260 template< class NI >
261 static constexpr std::size_t numIVAtVertex(const NI& nodalIndexSet)
262 { return 1; }
263
266 template< class IvIndexSetContainer,
267 class ScvfIndexMap,
268 class NodalIndexSet,
269 class FlipScvfIndexSet >
270 static void addIVIndexSets(IvIndexSetContainer& ivIndexSetContainer,
271 ScvfIndexMap& scvfIndexMap,
272 const NodalIndexSet& nodalIndexSet,
273 const FlipScvfIndexSet& flipScvfIndexSet)
274 {
275 // reuse the standard o-method's implementation of this
277 scvfIndexMap,
278 nodalIndexSet,
279 flipScvfIndexSet);
280 }
281
282private:
283 // pointer to cell stencil (in iv index set)
284 const Stencil * stencil_ = nullptr;
285
286 // Variables defining the local scope
287 std::array<Element, numScv> elements_;
288 std::array<LocalScvType, numScv> scvs_;
289 std::array<LocalScvfType, numScvf> scvfs_;
290 std::array<LocalFaceData, numScvf*2> localFaceData_;
291
292 // Dummy dirichlet data container (compatibility with dynamic o-iv)
293 std::array<DirichletData, 0> dirichletData_;
294};
295
296} // end namespace Dumux
297
298#endif
Define some often used mathematical functions.
Class for the index set within an interaction volume of the mpfa-o scheme.
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.
Class for the index sets of the dual grid in mpfa schemes.
MpfaMethods
The available mpfa schemes in Dumux.
Definition: methods.hh:34
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:68
General implementation of a data structure holding interaction volume-local information for a grid su...
Definition: localfacedata.hh:42
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:279
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 o-method's static interaction volume.
Definition: staticinteractionvolume.hh:126
static constexpr std::size_t numFaces()
returns the number of primary scvfs of this interaction volume
Definition: staticinteractionvolume.hh:214
void bind(const IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Sets up the local scope for a given iv index set.
Definition: staticinteractionvolume.hh:159
static constexpr MpfaMethods MpfaMethod
publicly state the mpfa-scheme this interaction volume is associated with
Definition: staticinteractionvolume.hh:155
static constexpr std::size_t numKnowns()
returns the number of (in this context) known solution values within this interaction volume
Definition: staticinteractionvolume.hh:222
static constexpr std::size_t numUnknowns()
returns the number of intermediate unknowns within this interaction volume
Definition: staticinteractionvolume.hh:218
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: staticinteractionvolume.hh:230
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the grid element corresponding to a given iv-local scv idx
Definition: staticinteractionvolume.hh:234
static constexpr std::size_t numIVAtVertex(const NI &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: staticinteractionvolume.hh:261
const std::array< LocalFaceData, numScvf *2 > & localFaceData() const
returns a reference to the container with the local face data
Definition: staticinteractionvolume.hh:246
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: staticinteractionvolume.hh:270
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: staticinteractionvolume.hh:256
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: staticinteractionvolume.hh:242
const std::array< DirichletData, 0 > & dirichletData() const
Definition: staticinteractionvolume.hh:251
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: staticinteractionvolume.hh:238
static constexpr std::size_t numScvs()
returns the number of scvs embedded in this interaction volume
Definition: staticinteractionvolume.hh:226
The default interaction volume traits class for the mpfa-o method with known size of the interaction ...
Definition: staticinteractionvolume.hh:65
static constexpr int numScvs
export the number of scvs in the interaction volumes
Definition: staticinteractionvolume.hh:104
typename NI::Traits::GridView GridView
export the type of grid view
Definition: staticinteractionvolume.hh:92
MVTraits MatVecTraits
export the matrix/vector traits to be used by the iv
Definition: staticinteractionvolume.hh:102
static constexpr int numScvfs
export the number of scvfs in the interaction volumes
Definition: staticinteractionvolume.hh:106
This does not work on surface grids.
Definition: staticinteractionvolume.hh:149
GridIndexType volVarIndex() const
Definition: staticinteractionvolume.hh:150
Classes for sub control entities of the mpfa-o method in the context of facet coupling.