3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
interactionvolumebase.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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_CC_MPFA_INTERACTIONVOLUMEBASE_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_INTERACTIONVOLUMEBASE_HH
26
27#include <vector>
28
29#include <dune/common/exceptions.hh>
30#include <dune/geometry/multilineargeometry.hh>
31
32namespace Dumux {
33
34/*
35 * \ingroup CCMpfaDiscretization
36 * \brief Type Traits to retrieve types associated with an implementation of a Dumux::CCMpfaInteractionVolume.
37 * You have to provide a traits class for every implementation of Dumux::CCMpfaInteractionVolume. Also,
38 * make sure to publicly export the traits class type in your interaction volume implementation.
39 * The traits should contain the following type definitions:
40 *
41 * \code
42 * //! export the type of grid view
43 * using GridView = ...;
44 * //! export the type used for local indices
45 * using IndexSet = ...;
46 * //! export the type of interaction-volume local scvs
47 * using LocalScvType = ...;
48 * //! export the type of interaction-volume local scvfs
49 * using LocalScvfType = ...;
50 * //! export the type of used for the iv-local face data
51 * using LocalFaceData = ...;
52 * //! export the matrix/vector type traits to be used by the iv
53 * using MatVecTraits = ...;
54 * //! export the type used for the assembly of the iv's local eq system
55 * using LocalAssembler = ...;
56 * \endcode
57 */
58
66template< class T >
68{
69 using GridView = typename T::GridView;
70 using Element = typename GridView::template Codim<0>::Entity;
71
72 using NodalStencilType = typename T::IndexSet::NodalGridStencilType;
73 using LocalIndexType = typename T::IndexSet::LocalIndexType;
74 using LocalScvType = typename T::LocalScvType;
75 using LocalScvfType = typename T::LocalScvfType;
76
77 using ScvGeometry = Dune::MultiLinearGeometry<typename LocalScvType::ctype,
78 LocalScvType::myDimension,
79 LocalScvType::worldDimension>;
80
81public:
83 using Traits = T;
84
86 template< class Problem, class FVElementGeometry >
87 void bind(const typename Traits::IndexSet& indexSet,
88 const Problem& problem,
89 const FVElementGeometry& fvGeometry)
90 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a bind() function"); }
91
93 std::size_t numFaces() const
94 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numFaces() function"); }
95
97 std::size_t numUnknowns() const
98 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numUnknowns() function"); }
99
101 std::size_t numKnowns() const
102 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numKnowns() function"); }
103
105 std::size_t numScvs() const
106 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numScvs() function"); }
107
109 template< class FVElementGeometry >
110 ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry& fvGeometry)
111 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a computeScvGeometry() function"); }
112
116 const std::vector<typename Traits::LocalFaceData>& localFaceData() const
117 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a localFaceData() function"); }
118
120 const NodalStencilType& stencil() const
121 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a stencil() function"); }
122
124 const LocalScvfType& localScvf(LocalIndexType ivLocalScvfIdx) const
125 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a localScvf() function"); }
126
128 const LocalScvType& localScv(LocalIndexType ivLocalScvIdx) const
129 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a localScv() function"); }
130
132 const Element& element(LocalIndexType ivLocalScvIdx) const
133 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide an element() function"); }
134
136 template< class NodalIndexSet >
137 static std::size_t numIVAtVertex(const NodalIndexSet& nodalIndexSet)
138 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numIVAtVertex() function"); }
139
142 template< class IvIndexSetContainer,
143 class ScvfIndexMap,
144 class NodalIndexSet,
145 class FlipScvfIndexSet >
146 static void addIVIndexSets(IvIndexSetContainer& ivIndexSetContainer,
147 ScvfIndexMap& scvfIndexMap,
148 const NodalIndexSet& nodalIndexSet,
149 const FlipScvfIndexSet& flipScvfIndexSet)
150 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide an addIVIndexSets() function"); }
151};
152
153} // end namespace Dumux
154
155#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:68
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: interactionvolumebase.hh:128
std::size_t numFaces() const
returns the number of "primary" scvfs of this interaction volume
Definition: interactionvolumebase.hh:93
void bind(const typename Traits::IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Prepares everything for the assembly.
Definition: interactionvolumebase.hh:87
T Traits
state the traits type publicly
Definition: interactionvolumebase.hh:83
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the element in which the scv with the given local idx is embedded in
Definition: interactionvolumebase.hh:132
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition: interactionvolumebase.hh:105
std::size_t numKnowns() const
returns the number of (in this context) known solution values within this interaction volume
Definition: interactionvolumebase.hh:101
static std::size_t numIVAtVertex(const NodalIndexSet &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: interactionvolumebase.hh:137
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: interactionvolumebase.hh:97
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: interactionvolumebase.hh:146
const std::vector< typename Traits::LocalFaceData > & localFaceData() const
Definition: interactionvolumebase.hh:116
const NodalStencilType & stencil() const
returns the cell-stencil of this interaction volume
Definition: interactionvolumebase.hh:120
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: interactionvolumebase.hh:124
ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: interactionvolumebase.hh:110