3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
staggered/gridfluxvariablescache.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_STAGGERED_GRID_FLUXVARSCACHE_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_GRID_FLUXVARSCACHE_HH
26
30
32
33namespace Dumux {
34
41template<class P, class FVC, class FVCF, int upwOrder>
43{
44 using Problem = P;
45 using FluxVariablesCache = FVC;
47
48 template<class GridFluxVariablesCache, bool cachingEnabled>
50 static constexpr int upwindSchemeOrder = upwOrder;
51};
52
57template<class Problem,
58 class FluxVariablesCache,
59 class FluxVariablesCacheFiller,
60 bool EnableGridFluxVariablesCache = false,
61 int upwindSchemeOrder = 1,
64
70template<class P, class FVC, class FVCF, int upwindSchemeOrder, class TheTraits>
71class StaggeredGridFluxVariablesCache<P, FVC, FVCF, true, upwindSchemeOrder, TheTraits>
72{
73 using Problem = typename TheTraits::Problem;
75
77 using FluxVariablesCacheFiller = typename TheTraits::FluxVariablesCacheFiller;
78public:
80 using Traits = TheTraits;
81
83 using FluxVariablesCache = typename Traits::FluxVariablesCache;
84 using Scalar = typename FluxVariablesCache::Scalar;
85
86 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
87
89 static constexpr bool cachingEnabled = true;
90
92 using LocalView = typename Traits::template LocalView<ThisType, cachingEnabled>;
93
94 StaggeredGridFluxVariablesCache(const Problem& problem, const std::string& paramGroup = "")
95 : problemPtr_(&problem)
96 , staggeredUpwindMethods_(paramGroup)
97 {}
98
99 // When global caching is enabled, precompute transmissibilities and stencils for all the scv faces
100 template<class GridGeometry, class GridVolumeVariables, class SolutionVector>
101 void update(const GridGeometry& gridGeometry,
102 const GridVolumeVariables& gridVolVars,
103 const SolutionVector& sol,
104 bool forceUpdate = false)
105 {
106 // only do the update if fluxes are solution dependent or if update is forced
107 if (FluxVariablesCacheFiller::isSolDependent || forceUpdate)
108 {
109 // instantiate helper class to fill the caches
110 // FluxVariablesCacheFiller filler(problem()); TODO: use proper ctor
111 FluxVariablesCacheFiller filler(problem());
112
113 fluxVarsCache_.resize(gridGeometry.numScvf());
114 for (const auto& element : elements(gridGeometry.gridView()))
115 {
116 // Prepare the geometries within the elements of the stencil
117 auto fvGeometry = localView(gridGeometry);
118 fvGeometry.bind(element);
119
120 auto elemVolVars = localView(gridVolVars);
121 elemVolVars.bind(element, fvGeometry, sol);
122
123 for (auto&& scvf : scvfs(fvGeometry))
124 {
125 filler.fill(*this, fluxVarsCache_[scvf.index()], element, fvGeometry, elemVolVars, scvf, forceUpdate);
126 }
127 }
128 }
129 }
130
133 {
134 return staggeredUpwindMethods_;
135 }
136
137 const Problem& problem() const
138 { return *problemPtr_; }
139
140 // access operators in the case of caching
141 const FluxVariablesCache& operator [](std::size_t scvfIdx) const
142 { return fluxVarsCache_[scvfIdx]; }
143
144 FluxVariablesCache& operator [](std::size_t scvfIdx)
145 { return fluxVarsCache_[scvfIdx]; }
146
147private:
148 const Problem* problemPtr_;
150
151 std::vector<FluxVariablesCache> fluxVarsCache_;
152 std::vector<std::size_t> globalScvfIndices_;
153};
154
160template<class P, class FVC, class FVCF, int upwindSchemeOrder, class TheTraits>
161class StaggeredGridFluxVariablesCache<P, FVC, FVCF, false, upwindSchemeOrder, TheTraits>
162{
163 using Problem = typename TheTraits::Problem;
165
167 using FluxVariablesCacheFiller = typename TheTraits::FluxVariablesCacheFiller;
168public:
170 using Traits = TheTraits;
171
173 using FluxVariablesCache = typename Traits::FluxVariablesCache;
174
176 using Scalar = typename FluxVariablesCache::Scalar;
177
179 static constexpr bool cachingEnabled = false;
180
182 using LocalView = typename Traits::template LocalView<ThisType, cachingEnabled>;
183
184 StaggeredGridFluxVariablesCache(const Problem& problem, const std::string& paramGroup = "")
185 : problemPtr_(&problem)
186 , staggeredUpwindMethods_(paramGroup)
187 {}
188
189 // When global caching is enabled, precompute transmissibilities and stencils for all the scv faces
190 template<class GridGeometry, class GridVolumeVariables, class SolutionVector>
191 void update(const GridGeometry& gridGeometry,
192 const GridVolumeVariables& gridVolVars,
193 const SolutionVector& sol,
194 bool forceUpdate = false) {}
195
196 const Problem& problem() const
197 { return *problemPtr_; }
198
201 {
202 return staggeredUpwindMethods_;
203 }
204
205private:
206 const Problem* problemPtr_;
208
209};
210
211} // end namespace Dumux
212
213#endif
Free function to get the local view of a grid cache object.
This file contains different higher order methods for approximating the velocity.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Base class for the stencil local flux variables cache for the staggered model.
Definition: staggered/elementfluxvariablescache.hh:41
Traits class to be used for the StaggeredGridVFluxVariablesCache.
Definition: staggered/gridfluxvariablescache.hh:43
static constexpr int upwindSchemeOrder
Definition: staggered/gridfluxvariablescache.hh:50
P Problem
Definition: staggered/gridfluxvariablescache.hh:44
FVCF FluxVariablesCacheFiller
Definition: staggered/gridfluxvariablescache.hh:46
FVC FluxVariablesCache
Definition: staggered/gridfluxvariablescache.hh:45
Flux variables cache class for staggered models.
Definition: staggered/gridfluxvariablescache.hh:63
Flux variables cache class for staggered models. Specialization in case of storing the flux cache.
Definition: staggered/gridfluxvariablescache.hh:72
typename FluxVariablesCache::Scalar Scalar
Definition: staggered/gridfluxvariablescache.hh:84
void update(const GridGeometry &gridGeometry, const GridVolumeVariables &gridVolVars, const SolutionVector &sol, bool forceUpdate=false)
Definition: staggered/gridfluxvariablescache.hh:101
TheTraits Traits
the flux var cache traits
Definition: staggered/gridfluxvariablescache.hh:80
typename Traits::template LocalView< ThisType, cachingEnabled > LocalView
export the type of the local view
Definition: staggered/gridfluxvariablescache.hh:92
StaggeredGridFluxVariablesCache(const Problem &problem, const std::string &paramGroup="")
Definition: staggered/gridfluxvariablescache.hh:94
typename Traits::FluxVariablesCache FluxVariablesCache
export the flux variable cache type
Definition: staggered/gridfluxvariablescache.hh:83
const StaggeredUpwindMethods< Scalar, upwindSchemeOrder > & staggeredUpwindMethods() const
Return the StaggeredUpwindMethods.
Definition: staggered/gridfluxvariablescache.hh:132
const Problem & problem() const
Definition: staggered/gridfluxvariablescache.hh:137
Flux variables cache class for staggered models. Specialization in case of not storing the flux cache...
Definition: staggered/gridfluxvariablescache.hh:162
StaggeredGridFluxVariablesCache(const Problem &problem, const std::string &paramGroup="")
Definition: staggered/gridfluxvariablescache.hh:184
typename Traits::FluxVariablesCache FluxVariablesCache
export the flux variable cache type
Definition: staggered/gridfluxvariablescache.hh:173
TheTraits Traits
the flux var cache traits
Definition: staggered/gridfluxvariablescache.hh:170
const Problem & problem() const
Definition: staggered/gridfluxvariablescache.hh:196
typename FluxVariablesCache::Scalar Scalar
the scalar type
Definition: staggered/gridfluxvariablescache.hh:176
const StaggeredUpwindMethods< Scalar, upwindSchemeOrder > & staggeredUpwindMethods() const
Return the UpwindingMethods.
Definition: staggered/gridfluxvariablescache.hh:200
typename Traits::template LocalView< ThisType, cachingEnabled > LocalView
export the type of the local view
Definition: staggered/gridfluxvariablescache.hh:182
void update(const GridGeometry &gridGeometry, const GridVolumeVariables &gridVolVars, const SolutionVector &sol, bool forceUpdate=false)
Definition: staggered/gridfluxvariablescache.hh:191
This file contains different higher order methods for approximating the velocity.
Definition: staggeredupwindmethods.hh:53