3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/fluxvariables.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_POROUSMEDIUMFLOW_FLUXVARIABLES_HH
26#define DUMUX_POROUSMEDIUMFLOW_FLUXVARIABLES_HH
27
28#include <bitset>
29#include <array>
30
34
35namespace Dumux {
36
46template<class TypeTag,
47 class UpScheme = UpwindScheme<GetPropType<TypeTag, Properties::GridGeometry>> >
49: public FluxVariablesBase<GetPropType<TypeTag, Properties::Problem>,
50 typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView,
51 typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView,
52 typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView>
53{
56
57 enum
58 {
59 numPhases = ModelTraits::numFluidPhases(),
60 numComponents = ModelTraits::numFluidComponents()
61 };
62
63public:
64 using UpwindScheme = UpScheme;
68
69 static constexpr bool enableAdvection = ModelTraits::enableAdvection();
70 static constexpr bool enableMolecularDiffusion = ModelTraits::enableMolecularDiffusion();
71 static constexpr bool enableEnergyBalance = ModelTraits::enableEnergyBalance();
72 static constexpr bool enableThermalNonEquilibrium = getPropValue<TypeTag, Properties::EnableThermalNonEquilibrium>();
73
76 {
77 advFluxIsCached_.reset();
78 advFluxBeforeUpwinding_.fill(0.0);
79 }
80
84 template<typename FunctionType>
85 Scalar advectiveFlux([[maybe_unused]] const int phaseIdx, [[maybe_unused]] const FunctionType& upwindTerm) const
86 {
87 if constexpr (enableAdvection)
88 {
89 if (!advFluxIsCached_[phaseIdx])
90 {
91
92 advFluxBeforeUpwinding_[phaseIdx] = AdvectionType::flux(this->problem(),
93 this->element(),
94 this->fvGeometry(),
95 this->elemVolVars(),
96 this->scvFace(),
97 phaseIdx,
98 this->elemFluxVarsCache());
99 advFluxIsCached_.set(phaseIdx, true);
100 }
101
103 return UpwindScheme::apply(*this, upwindTerm, advFluxBeforeUpwinding_[phaseIdx], phaseIdx);
104 }
105 else
106 return 0.0;
107 }
108
112 Dune::FieldVector<Scalar, numComponents> molecularDiffusionFlux([[maybe_unused]] const int phaseIdx) const
113 {
114 if constexpr (enableMolecularDiffusion)
115 return MolecularDiffusionType::flux(this->problem(),
116 this->element(),
117 this->fvGeometry(),
118 this->elemVolVars(),
119 this->scvFace(),
120 phaseIdx,
121 this->elemFluxVarsCache());
122 else
123 return Dune::FieldVector<Scalar, numComponents>(0.0);
124 }
125
130 Scalar heatConductionFlux() const
131 {
132 static_assert(!enableThermalNonEquilibrium, "This only works for thermal equilibrium");
133 if constexpr (enableEnergyBalance)
134 return HeatConductionType::flux(this->problem(),
135 this->element(),
136 this->fvGeometry(),
137 this->elemVolVars(),
138 this->scvFace(),
139 this->elemFluxVarsCache());
140 else
141 return 0.0;
142 }
143
148 Scalar heatConductionFlux([[maybe_unused]] const int phaseIdx) const
149 {
150 if constexpr (enableEnergyBalance)
151 return HeatConductionType::flux(this->problem(),
152 this->element(),
153 this->fvGeometry(),
154 this->elemVolVars(),
155 this->scvFace(),
156 phaseIdx,
157 this->elemFluxVarsCache());
158 else
159 return 0.0;
160 }
161
162private:
164 mutable std::bitset<numPhases> advFluxIsCached_;
165 mutable std::array<Scalar, numPhases> advFluxBeforeUpwinding_;
166};
167
168} // end namespace Dumux
169
170#endif
Base class for the flux variables living on a sub control volume face.
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Base class for the flux variables living on a sub control volume face.
Definition: fluxvariablesbase.hh:45
The porous medium flux variables class that computes advective / convective, molecular diffusive and ...
Definition: porousmediumflow/fluxvariables.hh:53
Scalar advectiveFlux(const int phaseIdx, const FunctionType &upwindTerm) const
Returns the advective flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:85
Dune::FieldVector< Scalar, numComponents > molecularDiffusionFlux(const int phaseIdx) const
Returns the diffusive fluxes computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:112
static constexpr bool enableAdvection
Definition: porousmediumflow/fluxvariables.hh:69
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: porousmediumflow/fluxvariables.hh:67
static constexpr bool enableThermalNonEquilibrium
Definition: porousmediumflow/fluxvariables.hh:72
static constexpr bool enableEnergyBalance
Definition: porousmediumflow/fluxvariables.hh:71
PorousMediumFluxVariables()
The constructor.
Definition: porousmediumflow/fluxvariables.hh:75
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition: porousmediumflow/fluxvariables.hh:66
UpScheme UpwindScheme
Definition: porousmediumflow/fluxvariables.hh:64
Scalar heatConductionFlux(const int phaseIdx) const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:148
Scalar heatConductionFlux() const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:130
GetPropType< TypeTag, Properties::AdvectionType > AdvectionType
Definition: porousmediumflow/fluxvariables.hh:65
static constexpr bool enableMolecularDiffusion
Definition: porousmediumflow/fluxvariables.hh:70
Declares all properties used in Dumux.
Base class for the upwind scheme.