3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
scalarfluxhelper.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_NAVIERSTOKES_SCALAR_BOUNDARY_FLUXHELPER_HH
25#define DUMUX_NAVIERSTOKES_SCALAR_BOUNDARY_FLUXHELPER_HH
26
27#include <dune/common/float_cmp.hh>
28#include <dune/common/std/type_traits.hh>
29#include <dumux/common/math.hh>
32
33namespace Dumux {
34
35#ifndef DOXYGEN
36namespace Detail {
37// helper structs and functions detecting if the VolumeVariables belong to a non-isothermal model
38template <class Indices>
39using NonisothermalDetector = decltype(std::declval<Indices>().energyEqIdx);
40
41template<class Indices>
42static constexpr bool isNonIsothermal()
43{ return Dune::Std::is_detected<NonisothermalDetector, Indices>::value; }
44
45} // end namespace Detail
46#endif
47
52template<class AdvectiveFlux>
54{
58 template<class VolumeVariables, class SubControlVolumeFace, class Scalar, class UpwindTerm>
59 static Scalar advectiveScalarUpwindFlux(const VolumeVariables& insideVolVars,
60 const VolumeVariables& outsideVolVars,
61 const SubControlVolumeFace& scvf,
62 const Scalar volumeFlux,
63 const Scalar upwindWeight,
64 UpwindTerm upwindTerm)
65 {
66 using std::signbit;
67 const bool insideIsUpstream = !signbit(volumeFlux);
68
69 const auto& upstreamVolVars = insideIsUpstream ? insideVolVars : outsideVolVars;
70 const auto& downstreamVolVars = insideIsUpstream ? outsideVolVars : insideVolVars;
71
72 return (upwindWeight * upwindTerm(upstreamVolVars) +
73 (1.0 - upwindWeight) * upwindTerm(downstreamVolVars))
74 * volumeFlux;
75 }
76
77 template<class Indices, class NumEqVector, class UpwindFunction>
79 const UpwindFunction& upwind)
80 {
81 // add advective fluxes based on physical type of model
82 AdvectiveFlux::addAdvectiveFlux(flux, upwind);
83
84 // for non-isothermal models, add the energy flux
85 if constexpr (Detail::isNonIsothermal<Indices>())
86 {
87 auto upwindTerm = [](const auto& volVars) { return volVars.density()*volVars.enthalpy(); };
88 flux[Indices::energyEqIdx] = upwind(upwindTerm);
89 }
90 }
91
96 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables>
97 static auto scalarOutflowFlux(const Problem& problem,
98 const Element& element,
99 const FVElementGeometry& fvGeometry,
100 const typename FVElementGeometry::SubControlVolumeFace& scvf,
101 const ElementVolumeVariables& elemVolVars,
102 typename ElementVolumeVariables::VolumeVariables::PrimaryVariables&& outsideBoundaryPriVars,
103 const typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type upwindWeight = 1.0)
104 {
105 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
106 using PrimaryVariables = typename VolumeVariables::PrimaryVariables;
107 using NumEqVector = PrimaryVariables;
108 NumEqVector flux;
109 const auto velocity = problem.faceVelocity(element,fvGeometry, scvf);
110 const auto volumeFlux = velocity * scvf.unitOuterNormal();
111 using std::signbit;
112 const bool insideIsUpstream = !signbit(volumeFlux);
113 const VolumeVariables& insideVolVars = elemVolVars[scvf.insideScvIdx()];
114 const VolumeVariables& outsideVolVars = [&]()
115 {
116 // only use the inside volVars for "true" outflow conditions (avoid constructing the outside volVars)
117 if (insideIsUpstream && Dune::FloatCmp::eq(upwindWeight, 1.0, 1e-6))
118 return insideVolVars;
119 else
120 {
121 // construct outside volVars from the given priVars for situations of flow reversal
122 VolumeVariables boundaryVolVars;
123 boundaryVolVars.update(elementSolution<FVElementGeometry>(std::forward<PrimaryVariables>(outsideBoundaryPriVars)),
124 problem,
125 element,
126 fvGeometry.scv(scvf.insideScvIdx()));
127 return boundaryVolVars;
128 }
129 }();
130
131 auto upwindFuntion = [&](const auto& upwindTerm)
132 {
133 return advectiveScalarUpwindFlux(insideVolVars, outsideVolVars, scvf, volumeFlux, upwindWeight, upwindTerm);
134 };
135
136 addModelSpecificAdvectiveFlux<typename VolumeVariables::Indices>(flux, upwindFuntion);
137
138 return flux;
139 }
140
146 template<class Problem, class Element, class FVElementGeometry, class ElementVolumeVariables>
147 static auto scalarOutflowFlux(const Problem& problem,
148 const Element& element,
149 const FVElementGeometry& fvGeometry,
150 const typename FVElementGeometry::SubControlVolumeFace& scvf,
151 const ElementVolumeVariables& elemVolVars)
152 {
153 using VolumeVariables = typename ElementVolumeVariables::VolumeVariables;
154 using NumEqVector = typename VolumeVariables::PrimaryVariables;
155 NumEqVector flux;
156 const auto velocity = problem.faceVelocity(element,fvGeometry, scvf);
157 const auto volumeFlux = velocity * scvf.unitOuterNormal();
158 using std::signbit;
159 const bool insideIsUpstream = !signbit(volumeFlux);
160 const VolumeVariables& insideVolVars = elemVolVars[scvf.insideScvIdx()];
161
162 if constexpr (VolumeVariables::FluidSystem::isCompressible(0/*phaseIdx*/) /*TODO viscosityIsConstant*/ || NumEqVector::size() > 1)
163 {
164 static const bool verbose = getParamFromGroup<bool>(problem.paramGroup(), "Flux.EnableOutflowReversalWarning", true);
165 using std::abs;
166 if (verbose && !insideIsUpstream && abs(volumeFlux) > 1e-10)
167 {
168 std::cout << "velo " << velocity << ", flux " << volumeFlux << std::endl;
169 std::cout << "\n ********** WARNING ********** \n\n"
170 "Outflow condition set at " << scvf.center() << " might be invalid due to flow reversal. "
171 "Consider using \n"
172 "outflowFlux(problem, element, fvGeometry, scvf, elemVolVars, outsideBoundaryPriVars, upwindWeight) \n"
173 "instead where you can specify primary variables for inflow situations.\n"
174 "\n ***************************** \n" << std::endl;
175 }
176 }
177
178 auto upwindFuntion = [&](const auto& upwindTerm)
179 {
180 return advectiveScalarUpwindFlux(insideVolVars, insideVolVars, scvf, volumeFlux, 1.0 /*upwindWeight*/, upwindTerm);
181 };
182
183 addModelSpecificAdvectiveFlux<typename VolumeVariables::Indices>(flux, upwindFuntion);
184
185 return flux;
186 }
187};
188
189} // end namespace Dumux
190
191#endif
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Element solution classes and factory functions.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
Definition: adapt.hh:29
Navier Stokes scalar boundary flux helper.
Definition: scalarfluxhelper.hh:54
static Scalar advectiveScalarUpwindFlux(const VolumeVariables &insideVolVars, const VolumeVariables &outsideVolVars, const SubControlVolumeFace &scvf, const Scalar volumeFlux, const Scalar upwindWeight, UpwindTerm upwindTerm)
Return the area-specific, weighted advective flux of a scalar quantity.
Definition: scalarfluxhelper.hh:59
static auto scalarOutflowFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars)
Return the area-specific outflow fluxes for all scalar balance equations. This should only be used of...
Definition: scalarfluxhelper.hh:147
static auto scalarOutflowFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, typename ElementVolumeVariables::VolumeVariables::PrimaryVariables &&outsideBoundaryPriVars, const typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type upwindWeight=1.0)
Return the area-specific outflow fluxes for all scalar balance equations. The values specified in out...
Definition: scalarfluxhelper.hh:97
static void addModelSpecificAdvectiveFlux(NumEqVector &flux, const UpwindFunction &upwind)
Definition: scalarfluxhelper.hh:78