24#ifndef DUMUX_NAVIERSTOKES_SCALAR_BOUNDARY_FLUXHELPER_HH
25#define DUMUX_NAVIERSTOKES_SCALAR_BOUNDARY_FLUXHELPER_HH
27#include <dune/common/float_cmp.hh>
28#include <dune/common/std/type_traits.hh>
38template <
class Indices>
39using NonisothermalDetector =
decltype(std::declval<Indices>().energyEqIdx);
41template<
class Indices>
42static constexpr bool isNonIsothermal()
43{
return Dune::Std::is_detected<NonisothermalDetector, Indices>::value; }
52template<
class AdvectiveFlux>
58 template<
class VolumeVariables,
class SubControlVolumeFace,
class Scalar,
class UpwindTerm>
60 const VolumeVariables& outsideVolVars,
61 const SubControlVolumeFace& scvf,
62 const Scalar volumeFlux,
63 const Scalar upwindWeight,
64 UpwindTerm upwindTerm)
67 const bool insideIsUpstream = !signbit(volumeFlux);
69 const auto& upstreamVolVars = insideIsUpstream ? insideVolVars : outsideVolVars;
70 const auto& downstreamVolVars = insideIsUpstream ? outsideVolVars : insideVolVars;
72 return (upwindWeight * upwindTerm(upstreamVolVars) +
73 (1.0 - upwindWeight) * upwindTerm(downstreamVolVars))
77 template<
class Indices,
class NumEqVector,
class UpwindFunction>
79 const UpwindFunction& upwind)
82 AdvectiveFlux::addAdvectiveFlux(flux, upwind);
85 if constexpr (Detail::isNonIsothermal<Indices>())
87 auto upwindTerm = [](
const auto& volVars) {
return volVars.density()*volVars.enthalpy(); };
88 flux[Indices::energyEqIdx] = upwind(upwindTerm);
96 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables>
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)
105 using VolumeVariables =
typename ElementVolumeVariables::VolumeVariables;
106 using PrimaryVariables =
typename VolumeVariables::PrimaryVariables;
109 const auto velocity = problem.faceVelocity(element,fvGeometry, scvf);
110 const auto volumeFlux = velocity * scvf.unitOuterNormal();
112 const bool insideIsUpstream = !signbit(volumeFlux);
113 const VolumeVariables& insideVolVars = elemVolVars[scvf.insideScvIdx()];
114 const VolumeVariables& outsideVolVars = [&]()
117 if (insideIsUpstream && Dune::FloatCmp::eq(upwindWeight, 1.0, 1e-6))
118 return insideVolVars;
122 VolumeVariables boundaryVolVars;
123 boundaryVolVars.update(elementSolution<FVElementGeometry>(std::forward<PrimaryVariables>(outsideBoundaryPriVars)),
126 fvGeometry.scv(scvf.insideScvIdx()));
127 return boundaryVolVars;
131 auto upwindFuntion = [&](
const auto& upwindTerm)
136 addModelSpecificAdvectiveFlux<typename VolumeVariables::Indices>(flux, upwindFuntion);
146 template<
class Problem,
class Element,
class FVElementGeometry,
class ElementVolumeVariables>
148 const Element& element,
149 const FVElementGeometry& fvGeometry,
150 const typename FVElementGeometry::SubControlVolumeFace& scvf,
151 const ElementVolumeVariables& elemVolVars)
153 using VolumeVariables =
typename ElementVolumeVariables::VolumeVariables;
154 using NumEqVector =
typename VolumeVariables::PrimaryVariables;
156 const auto velocity = problem.faceVelocity(element,fvGeometry, scvf);
157 const auto volumeFlux = velocity * scvf.unitOuterNormal();
159 const bool insideIsUpstream = !signbit(volumeFlux);
160 const VolumeVariables& insideVolVars = elemVolVars[scvf.insideScvIdx()];
162 if constexpr (VolumeVariables::FluidSystem::isCompressible(0) || NumEqVector::size() > 1)
164 static const bool verbose = getParamFromGroup<bool>(problem.paramGroup(),
"Flux.EnableOutflowReversalWarning",
true);
166 if (verbose && !insideIsUpstream && abs(volumeFlux) > 1e-10)
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. "
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;
178 auto upwindFuntion = [&](
const auto& upwindTerm)
183 addModelSpecificAdvectiveFlux<typename VolumeVariables::Indices>(flux, upwindFuntion);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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