version 3.10-dev
elementfacevariables.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_STAGGERED_ELEMENTFACEVARIABLES_HH
13#define DUMUX_DISCRETIZATION_STAGGERED_ELEMENTFACEVARIABLES_HH
14
15#include <algorithm>
16#include <cassert>
17#include <vector>
18#include <utility>
19
20namespace Dumux {
21
26template<class GridFaceVariables, bool cachingEnabled>
28{};
29
34template<class GFV>
35class StaggeredElementFaceVariables<GFV, /*cachingEnabled*/true>
36{
37public:
39 using GridFaceVariables = GFV;
40
42 using FaceVariables = typename GridFaceVariables::FaceVariables;
43
44 StaggeredElementFaceVariables(const GridFaceVariables& gridFaceVariables) : gridFaceVariablesPtr_(&gridFaceVariables) {}
45
47 template<class SubControlVolumeFace, typename std::enable_if_t<!std::is_integral<SubControlVolumeFace>::value, int> = 0>
48 const FaceVariables& operator [](const SubControlVolumeFace& scvf) const
49 { return gridFaceVariables().faceVars(scvf.index()); }
50
53 const FaceVariables& operator [](const std::size_t scvfIdx) const
54 { return gridFaceVariables().faceVars(scvfIdx); }
55
61 template<class FVElementGeometry, class SolutionVector>
62 StaggeredElementFaceVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
63 const FVElementGeometry& fvGeometry,
64 const SolutionVector& sol) &&
65 {
66 this->bind_(element, fvGeometry, sol);
67 return std::move(*this);
68 }
69
70 template<class FVElementGeometry, class SolutionVector>
71 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
72 const FVElementGeometry& fvGeometry,
73 const SolutionVector& sol) &
74 { this->bind_(element, fvGeometry, sol); }
75
81 template<class FVElementGeometry, class SolutionVector>
82 StaggeredElementFaceVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
83 const FVElementGeometry& fvGeometry,
84 const SolutionVector& sol) &&
85 {
86 this->bindElement_(element, fvGeometry, sol);
87 return std::move(*this);
88 }
89
90 template<class FVElementGeometry, class SolutionVector>
91 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
92 const FVElementGeometry& fvGeometry,
93 const SolutionVector& sol) &
94 { this->bindElement_(element, fvGeometry, sol); }
95
98 { return *gridFaceVariablesPtr_; }
99
100private:
101
104 template<class FVElementGeometry, class SolutionVector>
105 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
106 const FVElementGeometry& fvGeometry,
107 const SolutionVector& sol)
108 {}
109
112 template<class FVElementGeometry, class SolutionVector>
113 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
114 const FVElementGeometry& fvGeometry,
115 const SolutionVector& sol)
116 {}
117
118 const GridFaceVariables* gridFaceVariablesPtr_;
119};
120
125template<class GFV>
126class StaggeredElementFaceVariables<GFV, /*cachingEnabled*/false>
127{
128public:
130 using GridFaceVariables = GFV;
131
133 using FaceVariables = typename GridFaceVariables::FaceVariables;
134
135 StaggeredElementFaceVariables(const GridFaceVariables& globalFacesVars) : gridFaceVariablesPtr_(&globalFacesVars) {}
136
138 template<class SubControlVolumeFace, typename std::enable_if_t<!std::is_integral<SubControlVolumeFace>::value, int> = 0>
139 const FaceVariables& operator [](const SubControlVolumeFace& scvf) const
140 { return faceVariables_[scvf.localFaceIdx()]; }
141
143 const FaceVariables& operator [](const std::size_t scvfIdx) const
144 { return faceVariables_[getLocalIdx_(scvfIdx)]; }
145
147 template<class SubControlVolumeFace, typename std::enable_if_t<!std::is_integral<SubControlVolumeFace>::value, int> = 0>
148 FaceVariables& operator [](const SubControlVolumeFace& scvf)
149 { return faceVariables_[scvf.localFaceIdx()]; }
150
151 // operator for the access with an index
152 FaceVariables& operator [](const std::size_t scvfIdx)
153 { return faceVariables_[getLocalIdx_(scvfIdx)]; }
154
160 template<class FVElementGeometry, class SolutionVector>
161 StaggeredElementFaceVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
162 const FVElementGeometry& fvGeometry,
163 const SolutionVector& sol) &&
164 {
165 this->bind_(element, fvGeometry, sol);
166 return std::move(*this);
167 }
168
169 template<class FVElementGeometry, class SolutionVector>
170 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
171 const FVElementGeometry& fvGeometry,
172 const SolutionVector& sol) &
173 { this->bind_(element, fvGeometry, sol); }
174
180 template<class FVElementGeometry, class SolutionVector>
181 StaggeredElementFaceVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
182 const FVElementGeometry& fvGeometry,
183 const SolutionVector& sol) &&
184 {
185 this->bindElement_(element, fvGeometry, sol);
186 return std::move(*this);
187 }
188
189 template<class FVElementGeometry, class SolutionVector>
190 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
191 const FVElementGeometry& fvGeometry,
192 const SolutionVector& sol) &
193 { this->bindElement_(element, fvGeometry, sol); }
194
197 { return *gridFaceVariablesPtr_; }
198
199private:
200
203 template<class FVElementGeometry, class SolutionVector>
204 void bind_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
205 const FVElementGeometry& fvGeometry,
206 const SolutionVector& sol)
207 {
208 faceVariables_.resize(fvGeometry.numScvf());
209 faceVarIndices_.resize(fvGeometry.numScvf());
210
211 constexpr auto faceIdx = FVElementGeometry::GridGeometry::faceIdx();
212
213 for(auto&& scvf : scvfs(fvGeometry))
214 {
215 faceVariables_[scvf.localFaceIdx()].update(sol[faceIdx], gridFaceVariables().problem(), element, fvGeometry, scvf);
216 faceVarIndices_[scvf.localFaceIdx()] = scvf.index();
217 }
218 }
219
222 template<class FVElementGeometry, class SolutionVector>
223 void bindElement_(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
224 const FVElementGeometry& fvGeometry,
225 const SolutionVector& sol)
226 {
227 faceVariables_.resize(fvGeometry.numScvf());
228 faceVarIndices_.resize(fvGeometry.numScvf());
229
230 constexpr auto faceIdx = FVElementGeometry::GridGeometry::faceIdx();
231
232 for(auto&& scvf : scvfs(fvGeometry))
233 {
234 faceVariables_[scvf.localFaceIdx()].updateOwnFaceOnly(sol[faceIdx][scvf.dofIndex()]);
235 faceVarIndices_[scvf.localFaceIdx()] = scvf.index();
236 }
237 }
238
239 int getLocalIdx_(const int scvfIdx) const
240 {
241 auto it = std::find(faceVarIndices_.begin(), faceVarIndices_.end(), scvfIdx);
242 assert(it != faceVarIndices_.end() && "Could not find the current face variables for scvfIdx!");
243 return std::distance(faceVarIndices_.begin(), it);
244 }
245
246 const GridFaceVariables* gridFaceVariablesPtr_;
247 std::vector<std::size_t> faceVarIndices_;
248 std::vector<FaceVariables> faceVariables_;
249};
250
251} // end namespace
252
253#endif
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: elementfacevariables.hh:190
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: elementfacevariables.hh:170
StaggeredElementFaceVariables(const GridFaceVariables &globalFacesVars)
Definition: elementfacevariables.hh:135
typename GridFaceVariables::FaceVariables FaceVariables
export type of the volume variables
Definition: elementfacevariables.hh:133
StaggeredElementFaceVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: elementfacevariables.hh:181
const GridFaceVariables & gridFaceVariables() const
The global volume variables object we are a restriction of.
Definition: elementfacevariables.hh:196
GFV GridFaceVariables
export type of the grid volume variables
Definition: elementfacevariables.hh:130
StaggeredElementFaceVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: elementfacevariables.hh:161
StaggeredElementFaceVariables bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: elementfacevariables.hh:62
typename GridFaceVariables::FaceVariables FaceVariables
export type of the volume variables
Definition: elementfacevariables.hh:42
StaggeredElementFaceVariables(const GridFaceVariables &gridFaceVariables)
Definition: elementfacevariables.hh:44
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: elementfacevariables.hh:91
StaggeredElementFaceVariables bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: elementfacevariables.hh:82
const GridFaceVariables & gridFaceVariables() const
The global volume variables object we are a restriction of.
Definition: elementfacevariables.hh:97
GFV GridFaceVariables
export type of the grid volume variables
Definition: elementfacevariables.hh:39
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const SolutionVector &sol) &
Definition: elementfacevariables.hh:71
Base class for the face variables vector.
Definition: elementfacevariables.hh:28
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
Definition: adapt.hh:17