version 3.8
discretization/facecentered/staggered/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-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_FACECENTERED_STAGGERED_ELEMENT_FLUXVARSCACHE_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_FLUXVARSCACHE_HH
14
15#include <cstddef>
16#include <vector>
17#include <utility>
18
19namespace Dumux {
20
29template<class GFVC, bool cachingEnabled>
31
36template<class GFVC>
38{
39public:
42
44 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
45
47 : gridFluxVarsCachePtr_(&global) {}
48
54 template<class FVElementGeometry, class ElementVolumeVariables>
55 FaceCenteredStaggeredElementFluxVariablesCache bind(const typename FVElementGeometry::Element& element,
56 const FVElementGeometry& fvGeometry,
57 const ElementVolumeVariables& elemVolVars) &&
58 {
59 this->bind(element, fvGeometry, elemVolVars);
60 return std::move(*this);
61 }
62
63 // Function is called by the BoxLocalJacobian prior to flux calculations on the element.
64 // We assume the FVGeometries to be bound at this point
65 template<class FVElementGeometry, class ElementVolumeVariables>
66 void bind(const typename FVElementGeometry::Element& element,
67 const FVElementGeometry& fvGeometry,
68 const ElementVolumeVariables& elemVolVars) &
69 {
70 bindElement(element, fvGeometry, elemVolVars);
71 }
72
78 template<class FVElementGeometry, class ElementVolumeVariables>
79 FaceCenteredStaggeredElementFluxVariablesCache bindElement(const typename FVElementGeometry::Element& element,
80 const FVElementGeometry& fvGeometry,
81 const ElementVolumeVariables& elemVolVars) &&
82 {
83 this->bindElement(element, fvGeometry, elemVolVars);
84 return std::move(*this);
85 }
86
87 template<class FVElementGeometry, class ElementVolumeVariables>
88 void bindElement(const typename FVElementGeometry::Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars) &
91 {
92 eIdx_ = fvGeometry.gridGeometry().elementMapper().index(element);
93 }
94
100 template<class FVElementGeometry, class ElementVolumeVariables>
101 FaceCenteredStaggeredElementFluxVariablesCache bindScvf(const typename FVElementGeometry::Element& element,
102 const FVElementGeometry& fvGeometry,
103 const ElementVolumeVariables& elemVolVars,
104 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
105 {
106 this->bindScvf(element, fvGeometry, elemVolVars, scvf);
107 return std::move(*this);
108 }
109
110 template<class FVElementGeometry, class ElementVolumeVariables>
111 void bindScvf(const typename FVElementGeometry::Element& element,
112 const FVElementGeometry& fvGeometry,
113 const ElementVolumeVariables& elemVolVars,
114 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
115 {
116 bindElement(element, fvGeometry, elemVolVars);
117 }
118
119 // access operator
120 template<class SubControlVolumeFace>
121 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
122 { return gridFluxVarsCache()[scvf]; }
123
126 { return *gridFluxVarsCachePtr_; }
127
128private:
129 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
130 std::size_t eIdx_;
131};
132
137template<class GFVC>
139{
140public:
143
145 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
146
148 : gridFluxVarsCachePtr_(&global) {}
149
155 template<class FVElementGeometry, class ElementVolumeVariables>
156 FaceCenteredStaggeredElementFluxVariablesCache bind(const typename FVElementGeometry::Element& element,
157 const FVElementGeometry& fvGeometry,
158 const ElementVolumeVariables& elemVolVars) &&
159 {
160 this->bind(element, fvGeometry, elemVolVars);
161 return std::move(*this);
162 }
163
164 // Function is called by the BoxLocalJacobian prior to flux calculations on the element.
165 // We assume the FVGeometries to be bound at this point
166 template<class FVElementGeometry, class ElementVolumeVariables>
167 void bind(const typename FVElementGeometry::Element& element,
168 const FVElementGeometry& fvGeometry,
169 const ElementVolumeVariables& elemVolVars) &
170 {
171 bindElement(element, fvGeometry, elemVolVars);
172 }
173
179 template<class FVElementGeometry, class ElementVolumeVariables>
180 FaceCenteredStaggeredElementFluxVariablesCache bindElement(const typename FVElementGeometry::Element& element,
181 const FVElementGeometry& fvGeometry,
182 const ElementVolumeVariables& elemVolVars) &&
183 {
184 this->bindElement(element, fvGeometry, elemVolVars);
185 return std::move(*this);
186 }
187
188 template<class FVElementGeometry, class ElementVolumeVariables>
189 void bindElement(const typename FVElementGeometry::Element& element,
190 const FVElementGeometry& fvGeometry,
191 const ElementVolumeVariables& elemVolVars) &
192 {
193 // temporary resizing of the cache
194 fluxVarsCache_.resize(fvGeometry.numScvf());
195 // for (auto&& scvf : scvfs(fvGeometry))
196 // (*this)[scvf].update(gridFluxVarsCache().problem(), element, fvGeometry, elemVolVars, scvf); TODO
197 }
198
204 template<class FVElementGeometry, class ElementVolumeVariables>
205 FaceCenteredStaggeredElementFluxVariablesCache bindScvf(const typename FVElementGeometry::Element& element,
206 const FVElementGeometry& fvGeometry,
207 const ElementVolumeVariables& elemVolVars,
208 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
209 {
210 this->bindScvf(element, fvGeometry, elemVolVars, scvf);
211 return std::move(*this);
212 }
213
214 template<class FVElementGeometry, class ElementVolumeVariables>
215 void bindScvf(const typename FVElementGeometry::Element& element,
216 const FVElementGeometry& fvGeometry,
217 const ElementVolumeVariables& elemVolVars,
218 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
219 {
220 fluxVarsCache_.resize(fvGeometry.numScvf());
221 // (*this)[scvf].update(gridFluxVarsCache().problem(), element, fvGeometry, elemVolVars, scvf); TODO
222 }
223
224 // access operator
225 template<class SubControlVolumeFace>
226 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
227 { return fluxVarsCache_[scvf.index()]; }
228
229 // access operator
230 template<class SubControlVolumeFace>
231 FluxVariablesCache& operator [](const SubControlVolumeFace& scvf)
232 { return fluxVarsCache_[scvf.index()]; }
233
236 { return *gridFluxVarsCachePtr_; }
237
238private:
239 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
240 std::vector<FluxVariablesCache> fluxVarsCache_;
241};
242
243} // end namespace Dumux
244
245#endif
FaceCenteredStaggeredElementFluxVariablesCache bindScvf(const typename FVElementGeometry::Element &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/facecentered/staggered/elementfluxvariablescache.hh:205
FaceCenteredStaggeredElementFluxVariablesCache bind(const typename FVElementGeometry::Element &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/facecentered/staggered/elementfluxvariablescache.hh:156
void bindScvf(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:215
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:189
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:235
FaceCenteredStaggeredElementFluxVariablesCache(const GridFluxVariablesCache &global)
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:147
GFVC GridFluxVariablesCache
export the type of the grid flux variables cache
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:142
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:167
typename GFVC::FluxVariablesCache FluxVariablesCache
export the type of the flux variables cache
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:145
FaceCenteredStaggeredElementFluxVariablesCache bindElement(const typename FVElementGeometry::Element &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/facecentered/staggered/elementfluxvariablescache.hh:180
GFVC GridFluxVariablesCache
export the type of the grid flux variables cache
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:41
FaceCenteredStaggeredElementFluxVariablesCache(const GridFluxVariablesCache &global)
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:46
typename GFVC::FluxVariablesCache FluxVariablesCache
export the type of the flux variables cache
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:44
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:88
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:125
FaceCenteredStaggeredElementFluxVariablesCache bind(const typename FVElementGeometry::Element &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/facecentered/staggered/elementfluxvariablescache.hh:55
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:66
FaceCenteredStaggeredElementFluxVariablesCache bindElement(const typename FVElementGeometry::Element &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/facecentered/staggered/elementfluxvariablescache.hh:79
FaceCenteredStaggeredElementFluxVariablesCache bindScvf(const typename FVElementGeometry::Element &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/facecentered/staggered/elementfluxvariablescache.hh:101
void bindScvf(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:111
The flux variables caches for an element.
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:30
Definition: adapt.hh:17