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