3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_FACECENTERED_STAGGERED_ELEMENT_FLUXVARSCACHE_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_FLUXVARSCACHE_HH
26
27#include <cstddef>
28#include <vector>
29#include <utility>
30
31namespace Dumux {
32
41template<class GFVC, bool cachingEnabled>
43
48template<class GFVC>
50{
51public:
54
56 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
57
59 : gridFluxVarsCachePtr_(&global) {}
60
66 template<class FVElementGeometry, class ElementVolumeVariables>
67 FaceCenteredStaggeredElementFluxVariablesCache bind(const typename FVElementGeometry::Element& element,
68 const FVElementGeometry& fvGeometry,
69 const ElementVolumeVariables& elemVolVars) &&
70 {
71 this->bind(element, fvGeometry, elemVolVars);
72 return std::move(*this);
73 }
74
75 // Function is called by the BoxLocalJacobian prior to flux calculations on the element.
76 // We assume the FVGeometries to be bound at this point
77 template<class FVElementGeometry, class ElementVolumeVariables>
78 void bind(const typename FVElementGeometry::Element& element,
79 const FVElementGeometry& fvGeometry,
80 const ElementVolumeVariables& elemVolVars) &
81 {
82 bindElement(element, fvGeometry, elemVolVars);
83 }
84
90 template<class FVElementGeometry, class ElementVolumeVariables>
91 FaceCenteredStaggeredElementFluxVariablesCache bindElement(const typename FVElementGeometry::Element& 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 bindElement(const typename FVElementGeometry::Element& element,
101 const FVElementGeometry& fvGeometry,
102 const ElementVolumeVariables& elemVolVars) &
103 {
104 eIdx_ = fvGeometry.gridGeometry().elementMapper().index(element);
105 }
106
112 template<class FVElementGeometry, class ElementVolumeVariables>
113 FaceCenteredStaggeredElementFluxVariablesCache bindScvf(const typename FVElementGeometry::Element& element,
114 const FVElementGeometry& fvGeometry,
115 const ElementVolumeVariables& elemVolVars,
116 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
117 {
118 this->bindScvf(element, fvGeometry, elemVolVars, scvf);
119 return std::move(*this);
120 }
121
122 template<class FVElementGeometry, class ElementVolumeVariables>
123 void bindScvf(const typename FVElementGeometry::Element& element,
124 const FVElementGeometry& fvGeometry,
125 const ElementVolumeVariables& elemVolVars,
126 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
127 {
128 bindElement(element, fvGeometry, elemVolVars);
129 }
130
131 // access operator
132 template<class SubControlVolumeFace>
133 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
134 { return gridFluxVarsCache()[scvf]; }
135
138 { return *gridFluxVarsCachePtr_; }
139
140private:
141 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
142 std::size_t eIdx_;
143};
144
149template<class GFVC>
151{
152public:
155
157 using FluxVariablesCache = typename GFVC::FluxVariablesCache;
158
160 : gridFluxVarsCachePtr_(&global) {}
161
167 template<class FVElementGeometry, class ElementVolumeVariables>
168 FaceCenteredStaggeredElementFluxVariablesCache bind(const typename FVElementGeometry::Element& element,
169 const FVElementGeometry& fvGeometry,
170 const ElementVolumeVariables& elemVolVars) &&
171 {
172 this->bind(element, fvGeometry, elemVolVars);
173 return std::move(*this);
174 }
175
176 // Function is called by the BoxLocalJacobian prior to flux calculations on the element.
177 // We assume the FVGeometries to be bound at this point
178 template<class FVElementGeometry, class ElementVolumeVariables>
179 void bind(const typename FVElementGeometry::Element& element,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars) &
182 {
183 bindElement(element, fvGeometry, elemVolVars);
184 }
185
191 template<class FVElementGeometry, class ElementVolumeVariables>
192 FaceCenteredStaggeredElementFluxVariablesCache bindElement(const typename FVElementGeometry::Element& element,
193 const FVElementGeometry& fvGeometry,
194 const ElementVolumeVariables& elemVolVars) &&
195 {
196 this->bindElement(element, fvGeometry, elemVolVars);
197 return std::move(*this);
198 }
199
200 template<class FVElementGeometry, class ElementVolumeVariables>
201 void bindElement(const typename FVElementGeometry::Element& element,
202 const FVElementGeometry& fvGeometry,
203 const ElementVolumeVariables& elemVolVars) &
204 {
205 // temporary resizing of the cache
206 fluxVarsCache_.resize(fvGeometry.numScvf());
207 // for (auto&& scvf : scvfs(fvGeometry))
208 // (*this)[scvf].update(gridFluxVarsCache().problem(), element, fvGeometry, elemVolVars, scvf); TODO
209 }
210
216 template<class FVElementGeometry, class ElementVolumeVariables>
217 FaceCenteredStaggeredElementFluxVariablesCache bindScvf(const typename FVElementGeometry::Element& element,
218 const FVElementGeometry& fvGeometry,
219 const ElementVolumeVariables& elemVolVars,
220 const typename FVElementGeometry::SubControlVolumeFace& scvf) &&
221 {
222 this->bindScvf(element, fvGeometry, elemVolVars, scvf);
223 return std::move(*this);
224 }
225
226 template<class FVElementGeometry, class ElementVolumeVariables>
227 void bindScvf(const typename FVElementGeometry::Element& element,
228 const FVElementGeometry& fvGeometry,
229 const ElementVolumeVariables& elemVolVars,
230 const typename FVElementGeometry::SubControlVolumeFace& scvf) &
231 {
232 fluxVarsCache_.resize(fvGeometry.numScvf());
233 // (*this)[scvf].update(gridFluxVarsCache().problem(), element, fvGeometry, elemVolVars, scvf); TODO
234 }
235
236 // access operator
237 template<class SubControlVolumeFace>
238 const FluxVariablesCache& operator [](const SubControlVolumeFace& scvf) const
239 { return fluxVarsCache_[scvf.index()]; }
240
241 // access operator
242 template<class SubControlVolumeFace>
243 FluxVariablesCache& operator [](const SubControlVolumeFace& scvf)
244 { return fluxVarsCache_[scvf.index()]; }
245
248 { return *gridFluxVarsCachePtr_; }
249
250private:
251 const GridFluxVariablesCache* gridFluxVarsCachePtr_;
252 std::vector<FluxVariablesCache> fluxVarsCache_;
253};
254
255} // end namespace Dumux
256
257#endif
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
The flux variables caches for an element.
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:42
GFVC GridFluxVariablesCache
export the type of the grid flux variables cache
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:53
FaceCenteredStaggeredElementFluxVariablesCache(const GridFluxVariablesCache &global)
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:58
typename GFVC::FluxVariablesCache FluxVariablesCache
export the type of the flux variables cache
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:56
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:100
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:137
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:67
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:78
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:91
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:113
void bindScvf(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:123
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:217
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:168
void bindScvf(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:227
void bindElement(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:201
const GridFluxVariablesCache & gridFluxVarsCache() const
The global object we are a restriction of.
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:247
FaceCenteredStaggeredElementFluxVariablesCache(const GridFluxVariablesCache &global)
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:159
GFVC GridFluxVariablesCache
export the type of the grid flux variables cache
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:154
void bind(const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars) &
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:179
typename GFVC::FluxVariablesCache FluxVariablesCache
export the type of the flux variables cache
Definition: discretization/facecentered/staggered/elementfluxvariablescache.hh:157
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:192