3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
scalarfluxvariables.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_CONSERVATION_MODEL_FLUXVARIABLES_HH
25#define DUMUX_NAVIERSTOKES_SCALAR_CONSERVATION_MODEL_FLUXVARIABLES_HH
26
27#include <dumux/common/math.hh>
33
34namespace Dumux {
35
40template<class Problem,
41 class ModelTraits,
42 class FluxTypes,
43 class ElementVolumeVariables,
44 class ElementFluxVariablesCache,
45 class UpwindScheme = UpwindScheme<typename ProblemTraits<Problem>::GridGeometry>>
47: public FluxVariablesBase<Problem,
48 typename ProblemTraits<Problem>::GridGeometry::LocalView,
49 ElementVolumeVariables,
50 ElementFluxVariablesCache>
51{
52 using Scalar = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
53 using GridGeometry = typename ProblemTraits<Problem>::GridGeometry;
54 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
55 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
56 using FVElementGeometry = typename GridGeometry::LocalView;
58
59public:
60
62 {
63 static Scalar flux(const Problem& problem,
64 const Element& element,
65 const FVElementGeometry& fvGeometry,
66 const ElementVolumeVariables& elemVolVars,
67 const SubControlVolumeFace& scvf,
68 const int phaseIdx,
69 const ElementFluxVariablesCache& elemFluxVarsCache)
70 {
71 const auto velocity = problem.faceVelocity(element, fvGeometry, scvf);
72 const Scalar volumeFlux = velocity*scvf.unitOuterNormal()*Extrusion::area(fvGeometry, scvf)*extrusionFactor_(elemVolVars, scvf);
73 return volumeFlux;
74 }
75 };
76
80 template<typename FunctionType>
81 Scalar getAdvectiveFlux(const FunctionType& upwindTerm) const
82 {
83 if constexpr (ModelTraits::enableAdvection())
84 {
85 const auto& scvf = this->scvFace();
86 const auto velocity = this->problem().faceVelocity(this->element(), this->fvGeometry(), scvf);
87 const Scalar volumeFlux = velocity*scvf.unitOuterNormal()*Extrusion::area(this->fvGeometry(), scvf)*extrusionFactor_(this->elemVolVars(), scvf);
88 return UpwindScheme::apply(*this, upwindTerm, volumeFlux, 0/*phaseIdx*/);
89 }
90 else
91 return 0.0;
92 }
93
97 Scalar heatConductionFlux() const
98 {
99 if constexpr (ModelTraits::enableEnergyBalance())
100 {
101 using HeatConductionType = typename FluxTypes::HeatConductionType;
102 return HeatConductionType::flux(this->problem(),
103 this->element(),
104 this->fvGeometry(),
105 this->elemVolVars(),
106 this->scvFace(),
107 this->elemFluxVarsCache());
108 }
109 else
110 return 0.0;
111 }
112
116 Scalar heatAdvectionFlux() const
117 {
118 if constexpr (ModelTraits::enableEnergyBalance())
119 {
120 const auto upwindTerm = [](const auto& volVars) { return volVars.density() * volVars.enthalpy(); };
121 return getAdvectiveFlux(upwindTerm);
122 }
123 else
124 return 0.0;
125 }
126
130 Scalar heatFlux() const
131 {
133 }
134
138 template<class NumEqVector>
139 void addHeatFlux(NumEqVector& flux) const
140 {
141 if constexpr (ModelTraits::enableEnergyBalance())
142 flux[ModelTraits::Indices::energyEqIdx] = heatFlux();
143 }
144
145private:
146
147 static Scalar extrusionFactor_(const ElementVolumeVariables& elemVolVars, const SubControlVolumeFace& scvf)
148 {
149 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
150 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
151 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
152 }
153};
154
155} // end namespace Dumux
156
157#endif
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Base class for the flux variables living on a sub control volume face.
UpwindSchemeImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > UpwindScheme
The upwind scheme used for the advective fluxes. This depends on the chosen discretization method.
Definition: flux/upwindscheme.hh:42
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:69
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
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
std::decay_t< decltype(std::declval< Problem >().gridGeometry())> GridGeometry
Definition: common/typetraits/problem.hh:45
Base class for the flux variables living on a sub control volume face.
Definition: fluxvariablesbase.hh:45
const ProblemTraits< Problem >::GridGeometry::LocalView & fvGeometry() const
Definition: fluxvariablesbase.hh:78
The flux variables base class for scalar quantities balanced in the Navier-Stokes model.
Definition: scalarfluxvariables.hh:51
Scalar heatConductionFlux() const
Returns the conductive energy flux computed by the respective law.
Definition: scalarfluxvariables.hh:97
Scalar heatFlux() const
Returns the total energy flux.
Definition: scalarfluxvariables.hh:130
void addHeatFlux(NumEqVector &flux) const
Adds the energy flux to a given flux vector.
Definition: scalarfluxvariables.hh:139
Scalar heatAdvectionFlux() const
Returns the advective energy flux.
Definition: scalarfluxvariables.hh:116
Scalar getAdvectiveFlux(const FunctionType &upwindTerm) const
Returns the advective flux computed by the respective law.
Definition: scalarfluxvariables.hh:81
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: scalarfluxvariables.hh:63
Type traits for problem classes.
Base class for the upwind scheme.