version 3.11-dev
discretization/cvfe/hybrid/elementfluxvariablescache.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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_HYBRID_CVFE_ELEMENT_FLUXVARSCACHE_HH
13#define DUMUX_DISCRETIZATION_HYBRID_CVFE_ELEMENT_FLUXVARSCACHE_HH
14
15#include <array>
16#include <cassert>
17#include <cstddef>
18#include <optional>
19#include <ranges>
20#include <type_traits>
21#include <vector>
22#include <utility>
23
24#include <dumux/common/concepts/ipdata_.hh>
26
27namespace Dumux {
28
37template<class GFVC, bool cachingEnabled>
39
44template<class GFVC>
46{
47public:
50
52 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
53
55 : gridFluxVarsCachePtr_(&global) {}
56
57 // Function is called prior to flux calculations on the element.
58 // We assume the FVGeometries to be bound at this point
59 template<class FVElementGeometry, class ElementVolumeVariables>
60 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
61 const FVElementGeometry& fvGeometry,
62 const ElementVolumeVariables& elemVolVars) &
63 { bindElement(element, fvGeometry, elemVolVars); }
64
70 template<class FVElementGeometry, class ElementVolumeVariables>
71 HybridCVFEElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
72 const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars) &&
74 {
75 this->bind(element, fvGeometry, elemVolVars);
76 return std::move(*this);
77 }
78
79 template<class FVElementGeometry, class ElementVolumeVariables>
80 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
81 const FVElementGeometry& fvGeometry,
82 const ElementVolumeVariables& elemVolVars) &
83 { eIdx_ = fvGeometry.gridGeometry().elementMapper().index(element); }
84
90 template<class FVElementGeometry, class ElementVolumeVariables>
91 HybridCVFEElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
92 const FVElementGeometry& fvGeometry,
93 const ElementVolumeVariables& elemVolVars) &&
94 {
95 this->bindElement(element, fvGeometry, elemVolVars);
96 return std::move(*this);
97 }
98
99 template<class FVElementGeometry, class ElementVolumeVariables>
100 void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
101 const FVElementGeometry& fvGeometry,
102 const ElementVolumeVariables& elemVolVars,
103 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
104 { bindElement(element, fvGeometry, elemVolVars); }
105
111 template<class FVElementGeometry, class ElementVolumeVariables>
112 HybridCVFEElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
113 const FVElementGeometry& fvGeometry,
114 const ElementVolumeVariables& elemVolVars,
115 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
116 {
117 this->bindScvf(element, fvGeometry, elemVolVars, scvf);
118 return std::move(*this);
119 }
120
122 template<class FVElementGeometry, class ElementVolumeVariables>
123 void update(const typename FVElementGeometry::Element& element,
124 const FVElementGeometry& fvGeometry,
125 const ElementVolumeVariables& elemVolVars) {}
126
127 // access cache for given interpolation point data of an scvf quadrature point
128 template<Concept::ScvfQpIpData IpData>
129 const FluxVariablesCache& operator [](const IpData& ipData) const
130 {
131 return gridFluxVarsCache().cache(eIdx_, ipData.scvfIndex(), ipData.qpIndex());
132 }
133
134 // access cache for given interpolation point data of a boundary intersection quadrature point
135 template<Concept::IntersectionQpIpData IpData>
136 const FluxVariablesCache& operator [](const IpData& ipData) const
137 {
138 return gridFluxVarsCache().boundaryIntersectionCache(eIdx_, ipData.intersectionIndex(), ipData.qpIndex());
139 }
140
141 // access cache for given interpolation point data of an element quadrature point
142 template<Concept::QIpData IpData>
143 const FluxVariablesCache& operator [](const IpData& ipData) const
144 {
145 return gridFluxVarsCache().elementCache(eIdx_, ipData.qpIndex());
146 }
147
150 { return *gridFluxVarsCachePtr_; }
151
152private:
153 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
154 std::size_t eIdx_;
155};
156
161template<class GFVC>
163{
164 using GridGeometry = std::decay_t<decltype(std::declval<GFVC>().problem().gridGeometry())>;
165 using GridView = typename GridGeometry::GridView;
166
168 static constexpr std::size_t maxNumIntersections = GridView::dimension << 1;
169
170public:
173
175 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
176
178 : gridFluxVarsCachePtr_(&global) {}
179
180 // Function is called prior to flux calculations on the element
181 // We assume the FVGeometries to be bound at this point
182 template<class FVElementGeometry, class ElementVolumeVariables>
183 void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
184 const FVElementGeometry& fvGeometry,
185 const ElementVolumeVariables& elemVolVars) &
186 {
187 bindElement(element, fvGeometry, elemVolVars);
188 }
189
195 template<class FVElementGeometry, class ElementVolumeVariables>
196 HybridCVFEElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
197 const FVElementGeometry& fvGeometry,
198 const ElementVolumeVariables& elemVolVars) &&
199 {
200 this->bind(element, fvGeometry, elemVolVars);
201 return std::move(*this);
202 }
203
204 template<class FVElementGeometry, class ElementVolumeVariables>
205 void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
206 const FVElementGeometry& fvGeometry,
207 const ElementVolumeVariables& elemVolVars) &
208 {
209 // temporary resizing of the cache
210 fluxVarsCache_.resize(fvGeometry.numScvf());
211 for (auto&& scvf : scvfs(fvGeometry))
212 {
213 const auto quadRule = CVFE::quadratureRule(fvGeometry, scvf);
214 fluxVarsCache_[scvf.index()].resize(std::ranges::size(quadRule));
215 for (const auto& qpData : quadRule)
216 fluxVarsCache_[scvf.index()][qpData.ipData().qpIndex()].update(
217 gridFluxVarsCache().problem(), element, fvGeometry, elemVolVars, qpData.ipData()
218 );
219 }
220
221 // Resize element cache based on element quadrature rule
222 const auto elemQuadRule = CVFE::quadratureRule(fvGeometry, element);
223 elementCache_.resize(std::ranges::size(elemQuadRule));
224 for (const auto& qpData : elemQuadRule)
225 elementCache_[qpData.ipData().qpIndex()].update(gridFluxVarsCache().problem(), element, fvGeometry, elemVolVars, qpData.ipData());
226
227 // Set boundary intersection cache entries and reset non-boundary ones
228 for (const auto& intersection : intersections(fvGeometry.gridGeometry().gridView(), element))
229 {
230 const auto iIdx = intersection.indexInInside();
231 if (intersection.boundary())
232 {
233 auto& boundaryCache = boundaryIntersectionCache_[iIdx];
234 boundaryCache.emplace();
235 const auto quadRule = CVFE::quadratureRule(fvGeometry, intersection);
236 boundaryCache->resize(quadRule.size());
237 for (const auto& qpData : quadRule)
238 (*boundaryCache)[qpData.ipData().qpIndex()].update(gridFluxVarsCache().problem(),
239 element,
240 fvGeometry,
241 elemVolVars,
242 qpData.ipData());
243 }
244 else
245 boundaryIntersectionCache_[iIdx].reset();
246 }
247 }
248
254 template<class FVElementGeometry, class ElementVolumeVariables>
255 HybridCVFEElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
256 const FVElementGeometry& fvGeometry,
257 const ElementVolumeVariables& elemVolVars) &&
258 {
259 this->bindElement(element, fvGeometry, elemVolVars);
260 return std::move(*this);
261 }
262
263 template<class FVElementGeometry, class ElementVolumeVariables>
264 void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
265 const FVElementGeometry& fvGeometry,
266 const ElementVolumeVariables& elemVolVars,
267 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
268 {
269 fluxVarsCache_.resize(fvGeometry.numScvf());
270 const auto quadRule = CVFE::quadratureRule(fvGeometry, scvf);
271 fluxVarsCache_[scvf.index()].resize(std::ranges::size(quadRule));
272 for (const auto& qpData : quadRule)
273 fluxVarsCache_[scvf.index()][qpData.ipData().qpIndex()].update(
274 gridFluxVarsCache().problem(), element, fvGeometry, elemVolVars, qpData.ipData()
275 );
276 }
277
283 template<class FVElementGeometry, class ElementVolumeVariables>
284 HybridCVFEElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
285 const FVElementGeometry& fvGeometry,
286 const ElementVolumeVariables& elemVolVars,
287 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
288 {
289 this->bindScvf(element, fvGeometry, elemVolVars, scvf);
290 return std::move(*this);
291 }
292
297 template<class FVElementGeometry, class ElementVolumeVariables>
298 void update(const typename FVElementGeometry::Element& element,
299 const FVElementGeometry& fvGeometry,
300 const ElementVolumeVariables& elemVolVars)
301 {
302 if constexpr (FluxVariablesCache::isSolDependent)
303 bindElement(element, fvGeometry, elemVolVars);
304 }
305
306 // access cache for a given interpolation point data of an scvf quadrature point
307 template<Concept::ScvfQpIpData IpData>
308 const FluxVariablesCache& operator [](const IpData& ipData) const
309 { return fluxVarsCache_[ipData.scvfIndex()][ipData.qpIndex()]; }
310
311 // access cache for a given interpolation point data of an scvf quadrature point
312 template<Concept::ScvfQpIpData IpData>
313 FluxVariablesCache& operator [](const IpData& ipData)
314 { return fluxVarsCache_[ipData.scvfIndex()][ipData.qpIndex()]; }
315
316 // access cache for a given interpolation point data of an intersection quadrature point
317 template<Concept::IntersectionQpIpData IpData>
318 const FluxVariablesCache& operator [](const IpData& ipData) const
319 { return (*boundaryIntersectionCache_[ipData.intersectionIndex()])[ipData.qpIndex()]; }
320
321 // access cache for a given interpolation point data of an intersection quadrature point
322 template<Concept::IntersectionQpIpData IpData>
323 FluxVariablesCache& operator [](const IpData& ipData)
324 { return (*boundaryIntersectionCache_[ipData.intersectionIndex()])[ipData.qpIndex()]; }
325
326 // access cache for a given interpolation point data of an element quadrature point
327 template<Concept::QIpData IpData>
328 const FluxVariablesCache& operator [](const IpData& ipData) const
329 { return elementCache_[ipData.qpIndex()]; }
330
331 // access cache for a given interpolation point data of an element quadrature point
332 template<Concept::QIpData IpData>
333 FluxVariablesCache& operator [](const IpData& ipData)
334 { return elementCache_[ipData.qpIndex()]; }
335
338 { return *gridFluxVarsCachePtr_; }
339
340private:
341 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
342 std::vector<std::vector<FluxVariablesCache>> fluxVarsCache_;
343 std::vector<FluxVariablesCache> elementCache_;
344 std::array<std::optional<std::vector<FluxVariablesCache>>, maxNumIntersections> boundaryIntersectionCache_;
345};
346
347} // end namespace Dumux
348
349#endif
HybridCVFEElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:255
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:337
void update(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars)
Update the caches if the volume variables have changed and the cache is solution-dependent.
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:298
GFVC GridFluxVariablesCache
export the type of the grid flux variables cache
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:172
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:183
typename GFVC::FluxVariablesCache FluxVariablesCache
export the type of the flux variables cache
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:175
void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:264
HybridCVFEElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:196
HybridCVFEElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:284
HybridCVFEElementFluxVariablesCache(const GridFluxVariablesCache &global)
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:177
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:205
HybridCVFEElementFluxVariablesCache bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:112
typename GFVC::FluxVariablesCache FluxVariablesCache
export the type of the flux variables cache
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:52
GFVC GridFluxVariablesCache
export the type of the grid flux variables cache
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:49
void bindScvf(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:100
void bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:80
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:149
HybridCVFEElementFluxVariablesCache(const GridFluxVariablesCache &global)
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:54
void update(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars)
Specialization for the global caching being enabled - do nothing here.
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:123
HybridCVFEElementFluxVariablesCache bindElement(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:91
HybridCVFEElementFluxVariablesCache bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:71
void bind(const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:60
The flux variables caches for an element when using hybrid CVFE discretizations.
Definition: discretization/cvfe/hybrid/elementfluxvariablescache.hh:38
auto quadratureRule(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolume &scv, QuadratureRules::MidpointQuadrature)
Midpoint quadrature for scv.
Definition: quadraturerules.hh:148
Definition: adapt.hh:17
Quadrature rules over sub-control volumes and sub-control volume faces.