3.1-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 = getParam<Scalar>("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 if (std::signbit(flux)) // if sign of flux is negative
66 return flux*(upwindWeight*upwindTerm(outsideVolVars)
67 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
68 else
69 return flux*(upwindWeight*upwindTerm(insideVolVars)
70 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
71 }
72};
73
75template<class GridGeometry>
77{
78 using GridView = typename GridGeometry::GridView;
79 static constexpr int dim = GridView::dimension;
80 static constexpr int dimWorld = GridView::dimensionworld;
81
82public:
83 // For surface and network grids (dim < dimWorld) we have to do a special upwind scheme
84 template<class FluxVariables, class UpwindTermFunction, class Scalar, int d = dim, int dw = dimWorld>
85 static typename std::enable_if<(d < dw), Scalar>::type
86 apply(const FluxVariables& fluxVars,
87 const UpwindTermFunction& upwindTerm,
88 Scalar flux, int phaseIdx)
89 {
90 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
91
92 // the volume variables of the inside sub-control volume
93 const auto& scvf = fluxVars.scvFace();
94 const auto& elemVolVars = fluxVars.elemVolVars();
95 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
96
97 // check if this is a branching point
98 if (scvf.numOutsideScvs() > 1)
99 {
100 // more complicated upwind scheme
101 // we compute a flux-weighted average of all inflowing branches
102 Scalar branchingPointUpwindTerm = 0.0;
103 Scalar sumUpwindFluxes = 0.0;
104
105 // if the inside flux is positive (outflow) do fully upwind and return flux
106 if (!std::signbit(flux))
107 return upwindTerm(insideVolVars)*flux;
108 else
109 sumUpwindFluxes += flux;
110
111 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
112 {
113 // compute the outside flux
114 const auto& fvGeometry = fluxVars.fvGeometry();
115 const auto outsideScvIdx = scvf.outsideScvIdx(i);
116 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
117 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
118
119 using AdvectionType = typename FluxVariables::AdvectionType;
120 const auto outsideFlux = AdvectionType::flux(fluxVars.problem(),
121 outsideElement,
122 fvGeometry,
123 elemVolVars,
124 flippedScvf,
125 phaseIdx,
126 fluxVars.elemFluxVarsCache());
127
128 if (!std::signbit(outsideFlux))
129 branchingPointUpwindTerm += upwindTerm(elemVolVars[outsideScvIdx])*outsideFlux;
130 else
131 sumUpwindFluxes += outsideFlux;
132 }
133
134 // the flux might be zero
135 if (sumUpwindFluxes != 0.0)
136 branchingPointUpwindTerm /= -sumUpwindFluxes;
137 else
138 branchingPointUpwindTerm = 0.0;
139
140 // upwind scheme (always do fully upwind at branching points)
141 // a weighting here would lead to an error since the derivation is based on a fully upwind scheme
142 // TODO How to implement a weight of e.g. 0.5
143 if (std::signbit(flux))
144 return flux*branchingPointUpwindTerm;
145 else
146 return flux*upwindTerm(insideVolVars);
147 }
148 // non-branching points and boundaries
149 else
150 {
151 // upwind scheme
152 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
153 if (std::signbit(flux))
154 return flux*(upwindWeight*upwindTerm(outsideVolVars)
155 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
156 else
157 return flux*(upwindWeight*upwindTerm(insideVolVars)
158 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
159 }
160 }
161
162 // For grids with dim == dimWorld we use a simple upwinding scheme
163 template<class FluxVariables, class UpwindTermFunction, class Scalar, int d = dim, int dw = dimWorld>
164 static typename std::enable_if<(d == dw), Scalar>::type
165 apply(const FluxVariables& fluxVars,
166 const UpwindTermFunction& upwindTerm,
167 Scalar flux, int phaseIdx)
168 {
169 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
170
171 const auto& scvf = fluxVars.scvFace();
172 const auto& elemVolVars = fluxVars.elemVolVars();
173 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
174 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
175
176 if (std::signbit(flux)) // if sign of flux is negative
177 return flux*(upwindWeight*upwindTerm(outsideVolVars)
178 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
179 else
180 return flux*(upwindWeight*upwindTerm(insideVolVars)
181 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
182 }
183};
184
186template<class GridGeometry>
188: public UpwindSchemeImpl<GridGeometry, DiscretizationMethod::cctpfa> {};
189
190} // end namespace Dumux
191
192#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
make the local view function available whenever we use the grid geometry
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:86
static std::enable_if<(d==dw), Scalar >::type apply(const FluxVariables &fluxVars, const UpwindTermFunction &upwindTerm, Scalar flux, int phaseIdx)
Definition: flux/upwindscheme.hh:165