3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/upwindscheme.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_UPWINDSCHEME_HH
25#define DUMUX_DISCRETIZATION_UPWINDSCHEME_HH
26
29
30namespace Dumux {
31
33template<class GridGeometry, DiscretizationMethod discMethod>
35
41template<class GridGeometry>
43
45template<class GridGeometry>
47{
48public:
49 // applies a simple upwind scheme to the precalculated advective flux
50 template<class FluxVariables, class UpwindTermFunction, class Scalar>
51 static Scalar apply(const FluxVariables& fluxVars,
52 const UpwindTermFunction& upwindTerm,
53 Scalar flux, int phaseIdx)
54 {
55 // TODO: pass this from outside?
56 static const Scalar upwindWeight = getParamFromGroup<Scalar>(fluxVars.problem().paramGroup(), "Flux.UpwindWeight");
57
58 const auto& elemVolVars = fluxVars.elemVolVars();
59 const auto& scvf = fluxVars.scvFace();
60 const auto& insideScv = fluxVars.fvGeometry().scv(scvf.insideScvIdx());
61 const auto& outsideScv = fluxVars.fvGeometry().scv(scvf.outsideScvIdx());
62 const auto& insideVolVars = elemVolVars[insideScv];
63 const auto& outsideVolVars = elemVolVars[outsideScv];
64
65 using std::signbit;
66 if (signbit(flux)) // if sign of flux is negative
67 return flux*(upwindWeight*upwindTerm(outsideVolVars)
68 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
69 else
70 return flux*(upwindWeight*upwindTerm(insideVolVars)
71 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
72 }
73};
74
76template<class GridGeometry>
78{
79 using GridView = typename GridGeometry::GridView;
80 static constexpr int dim = GridView::dimension;
81 static constexpr int dimWorld = GridView::dimensionworld;
82
83public:
84 // For surface and network grids (dim < dimWorld) we have to do a special upwind scheme
85 template<class FluxVariables, class UpwindTermFunction, class Scalar, int d = dim, int dw = dimWorld>
86 static typename std::enable_if<(d < dw), Scalar>::type
87 apply(const FluxVariables& fluxVars,
88 const UpwindTermFunction& upwindTerm,
89 Scalar flux, int phaseIdx)
90 {
91 using std::signbit;
92 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
93
94 // the volume variables of the inside sub-control volume
95 const auto& scvf = fluxVars.scvFace();
96 const auto& elemVolVars = fluxVars.elemVolVars();
97 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
98
99 // check if this is a branching point
100 if (scvf.numOutsideScvs() > 1)
101 {
102 // more complicated upwind scheme
103 // we compute a flux-weighted average of all inflowing branches
104 Scalar branchingPointUpwindTerm = 0.0;
105 Scalar sumUpwindFluxes = 0.0;
106
107 // if the inside flux is positive (outflow) do fully upwind and return flux
108 if (!signbit(flux))
109 return upwindTerm(insideVolVars)*flux;
110 else
111 sumUpwindFluxes += flux;
112
113 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
114 {
115 // compute the outside flux
116 const auto& fvGeometry = fluxVars.fvGeometry();
117 const auto outsideScvIdx = scvf.outsideScvIdx(i);
118 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
119 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
120
121 using AdvectionType = typename FluxVariables::AdvectionType;
122 const auto outsideFlux = AdvectionType::flux(fluxVars.problem(),
123 outsideElement,
124 fvGeometry,
125 elemVolVars,
126 flippedScvf,
127 phaseIdx,
128 fluxVars.elemFluxVarsCache());
129
130 if (!signbit(outsideFlux))
131 branchingPointUpwindTerm += upwindTerm(elemVolVars[outsideScvIdx])*outsideFlux;
132 else
133 sumUpwindFluxes += outsideFlux;
134 }
135
136 // the flux might be zero
137 if (sumUpwindFluxes != 0.0)
138 branchingPointUpwindTerm /= -sumUpwindFluxes;
139 else
140 branchingPointUpwindTerm = 0.0;
141
142 // upwind scheme (always do fully upwind at branching points)
143 // a weighting here would lead to an error since the derivation is based on a fully upwind scheme
144 // TODO How to implement a weight of e.g. 0.5
145 if (signbit(flux))
146 return flux*branchingPointUpwindTerm;
147 else
148 return flux*upwindTerm(insideVolVars);
149 }
150 // non-branching points and boundaries
151 else
152 {
153 // upwind scheme
154 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
155 if (signbit(flux))
156 return flux*(upwindWeight*upwindTerm(outsideVolVars)
157 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
158 else
159 return flux*(upwindWeight*upwindTerm(insideVolVars)
160 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
161 }
162 }
163
164 // For grids with dim == dimWorld we use a simple upwinding scheme
165 template<class FluxVariables, class UpwindTermFunction, class Scalar, int d = dim, int dw = dimWorld>
166 static typename std::enable_if<(d == dw), Scalar>::type
167 apply(const FluxVariables& fluxVars,
168 const UpwindTermFunction& upwindTerm,
169 Scalar flux, int phaseIdx)
170 {
171 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
172
173 const auto& scvf = fluxVars.scvFace();
174 const auto& elemVolVars = fluxVars.elemVolVars();
175 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
176 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
177
178 using std::signbit;
179 if (signbit(flux)) // if sign of flux is negative
180 return flux*(upwindWeight*upwindTerm(outsideVolVars)
181 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
182 else
183 return flux*(upwindWeight*upwindTerm(insideVolVars)
184 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
185 }
186};
187
189template<class GridGeometry>
191: public UpwindSchemeImpl<GridGeometry, DiscretizationMethod::cctpfa> {};
192
193} // end namespace Dumux
194
195#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Definition: adapt.hh:29
Forward declaration of the upwind scheme implementation.
Definition: flux/upwindscheme.hh:34
static Scalar apply(const FluxVariables &fluxVars, const UpwindTermFunction &upwindTerm, Scalar flux, int phaseIdx)
Definition: flux/upwindscheme.hh:51
static std::enable_if<(d< dw), Scalar >::type apply(const FluxVariables &fluxVars, const UpwindTermFunction &upwindTerm, Scalar flux, int phaseIdx)
Definition: flux/upwindscheme.hh:87
static std::enable_if<(d==dw), Scalar >::type apply(const FluxVariables &fluxVars, const UpwindTermFunction &upwindTerm, Scalar flux, int phaseIdx)
Definition: flux/upwindscheme.hh:167