version 3.8
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// 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_DISCRETIZATION_CC_MPFA_INTERACTIONVOLUMEBASE_HH
13#define DUMUX_DISCRETIZATION_CC_MPFA_INTERACTIONVOLUMEBASE_HH
14
15#include <vector>
16
17#include <dune/common/exceptions.hh>
18#include <dune/geometry/multilineargeometry.hh>
19
20namespace Dumux {
21
22/*
23 * \ingroup CCMpfaDiscretization
24 * \brief Type Traits to retrieve types associated with an implementation of a Dumux::CCMpfaInteractionVolume.
25 * You have to provide a traits class for every implementation of Dumux::CCMpfaInteractionVolume. Also,
26 * make sure to publicly export the traits class type in your interaction volume implementation.
27 * The traits should contain the following type definitions:
28 *
29 * \code
30 * //! export the type of grid view
31 * using GridView = ...;
32 * //! export the type used for local indices
33 * using IndexSet = ...;
34 * //! export the type of interaction-volume local scvs
35 * using LocalScvType = ...;
36 * //! export the type of interaction-volume local scvfs
37 * using LocalScvfType = ...;
38 * //! export the type of used for the iv-local face data
39 * using LocalFaceData = ...;
40 * //! export the matrix/vector type traits to be used by the iv
41 * using MatVecTraits = ...;
42 * //! export the type used for the assembly of the iv's local eq system
43 * using LocalAssembler = ...;
44 * \endcode
45 */
46
54template< class T >
56{
57 using GridView = typename T::GridView;
58 using Element = typename GridView::template Codim<0>::Entity;
59
60 using NodalStencilType = typename T::IndexSet::NodalGridStencilType;
61 using LocalIndexType = typename T::IndexSet::LocalIndexType;
62 using LocalScvType = typename T::LocalScvType;
63 using LocalScvfType = typename T::LocalScvfType;
64
65 using ScvGeometry = Dune::MultiLinearGeometry<typename LocalScvType::ctype,
66 LocalScvType::myDimension,
67 LocalScvType::worldDimension>;
68
69public:
71 using Traits = T;
72
74 template< class Problem, class FVElementGeometry >
75 void bind(const typename Traits::IndexSet& indexSet,
76 const Problem& problem,
77 const FVElementGeometry& fvGeometry)
78 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a bind() function"); }
79
81 std::size_t numFaces() const
82 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numFaces() function"); }
83
85 std::size_t numUnknowns() const
86 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numUnknowns() function"); }
87
89 std::size_t numKnowns() const
90 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numKnowns() function"); }
91
93 std::size_t numScvs() const
94 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numScvs() function"); }
95
97 template< class FVElementGeometry >
98 ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry& fvGeometry)
99 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a computeScvGeometry() function"); }
100
104 const std::vector<typename Traits::LocalFaceData>& localFaceData() const
105 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a localFaceData() function"); }
106
108 const NodalStencilType& stencil() const
109 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a stencil() function"); }
110
112 const LocalScvfType& localScvf(LocalIndexType ivLocalScvfIdx) const
113 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a localScvf() function"); }
114
116 const LocalScvType& localScv(LocalIndexType ivLocalScvIdx) const
117 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a localScv() function"); }
118
120 const Element& element(LocalIndexType ivLocalScvIdx) const
121 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide an element() function"); }
122
124 template< class NodalIndexSet >
125 static std::size_t numIVAtVertex(const NodalIndexSet& nodalIndexSet)
126 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide a numIVAtVertex() function"); }
127
130 template< class IvIndexSetContainer,
131 class ScvfIndexMap,
132 class NodalIndexSet,
133 class FlipScvfIndexSet >
134 static void addIVIndexSets(IvIndexSetContainer& ivIndexSetContainer,
135 ScvfIndexMap& scvfIndexMap,
136 const NodalIndexSet& nodalIndexSet,
137 const FlipScvfIndexSet& flipScvfIndexSet)
138 { DUNE_THROW(Dune::NotImplemented, "Interaction volume does not provide an addIVIndexSets() function"); }
139};
140
141} // end namespace Dumux
142
143#endif
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:56
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: interactionvolumebase.hh:116
std::size_t numFaces() const
returns the number of "primary" scvfs of this interaction volume
Definition: interactionvolumebase.hh:81
void bind(const typename Traits::IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Prepares everything for the assembly.
Definition: interactionvolumebase.hh:75
T Traits
state the traits type publicly
Definition: interactionvolumebase.hh:71
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the element in which the scv with the given local idx is embedded in
Definition: interactionvolumebase.hh:120
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition: interactionvolumebase.hh:93
std::size_t numKnowns() const
returns the number of (in this context) known solution values within this interaction volume
Definition: interactionvolumebase.hh:89
static std::size_t numIVAtVertex(const NodalIndexSet &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: interactionvolumebase.hh:125
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: interactionvolumebase.hh:85
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: interactionvolumebase.hh:134
const std::vector< typename Traits::LocalFaceData > & localFaceData() const
Definition: interactionvolumebase.hh:104
const NodalStencilType & stencil() const
returns the cell-stencil of this interaction volume
Definition: interactionvolumebase.hh:108
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: interactionvolumebase.hh:112
ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: interactionvolumebase.hh:98
Definition: adapt.hh:17