3.1-git
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 UpwindScheme = 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{
57 using FVElementGeometry = typename GridGeometry::LocalView;
58 using GridView = typename GridGeometry::GridView;
59 using Element = typename GridView::template Codim<0>::Entity;
60 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
61
63 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
64 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
65
66 enum
67 {
68 numPhases = ModelTraits::numFluidPhases(),
69 numComponents = ModelTraits::numFluidComponents()
70 };
71
72public:
76
77 static constexpr bool enableAdvection = ModelTraits::enableAdvection();
78 static constexpr bool enableMolecularDiffusion = ModelTraits::enableMolecularDiffusion();
79 static constexpr bool enableEnergyBalance = ModelTraits::enableEnergyBalance();
80 static constexpr bool enableThermalNonEquilibrium = getPropValue<TypeTag, Properties::EnableThermalNonEquilibrium>();
81
84 {
85 advFluxIsCached_.reset();
86 advFluxBeforeUpwinding_.fill(0.0);
87 }
88
93 template<typename FunctionType, bool enable = enableAdvection, std::enable_if_t<enable, int> = 0>
94 Scalar advectiveFlux(const int phaseIdx, const FunctionType& upwindTerm) const
95 {
96 if (!advFluxIsCached_[phaseIdx])
97 {
98
99 advFluxBeforeUpwinding_[phaseIdx] = AdvectionType::flux(this->problem(),
100 this->element(),
101 this->fvGeometry(),
102 this->elemVolVars(),
103 this->scvFace(),
104 phaseIdx,
105 this->elemFluxVarsCache());
106 advFluxIsCached_.set(phaseIdx, true);
107 }
108
110 return UpwindScheme::apply(*this, upwindTerm, advFluxBeforeUpwinding_[phaseIdx], phaseIdx);
111 }
112
117 template<typename FunctionType, bool enable = enableAdvection, typename std::enable_if_t<!enable, int> = 0>
118 Scalar advectiveFlux(const int phaseIdx, const FunctionType& upwindTerm) const
119 {
120 return 0.0;
121 }
122
128 template<bool enable = enableMolecularDiffusion, typename std::enable_if_t<enable, int> = 0>
129 Dune::FieldVector<Scalar, numComponents> molecularDiffusionFlux(const int phaseIdx) const
130 {
131 return MolecularDiffusionType::flux(this->problem(),
132 this->element(),
133 this->fvGeometry(),
134 this->elemVolVars(),
135 this->scvFace(),
136 phaseIdx,
137 this->elemFluxVarsCache());
138 }
139
145 template<bool enable = enableMolecularDiffusion, typename std::enable_if_t<!enable, int> = 0>
146 Dune::FieldVector<Scalar, numComponents> molecularDiffusionFlux(const int phaseIdx) const
147 {
148 return Dune::FieldVector<Scalar, numComponents>(0.0);
149 }
150
156 template<bool enable = enableEnergyBalance && !enableThermalNonEquilibrium, typename std::enable_if_t<enable, int> = 0>
157 Scalar heatConductionFlux() const
158 {
159 return HeatConductionType::flux(this->problem(),
160 this->element(),
161 this->fvGeometry(),
162 this->elemVolVars(),
163 this->scvFace(),
164 this->elemFluxVarsCache());
165 }
166
172 template<bool enable = enableEnergyBalance && enableThermalNonEquilibrium, typename std::enable_if_t<enable, int> = 0>
173 Scalar heatConductionFlux(const int phaseIdx) const
174 {
175 return HeatConductionType::flux(this->problem(),
176 this->element(),
177 this->fvGeometry(),
178 this->elemVolVars(),
179 this->scvFace(),
180 phaseIdx,
181 this->elemFluxVarsCache());
182 }
183
189 template<bool enable = enableEnergyBalance, typename std::enable_if_t<!enable, int> = 0>
190 Scalar heatConductionFlux(const int phaseIdx) const
191 {
192 return 0.0;
193 }
194
195private:
197 mutable std::bitset<numPhases> advFluxIsCached_;
198 mutable std::array<Scalar, numPhases> advFluxBeforeUpwinding_;
199};
200
201} // end namespace Dumux
202
203#endif
Base class for the flux variables living on a sub control volume face.
UpwindSchemeImpl< GridGeometry, GridGeometry::discMethod > UpwindScheme
The upwind scheme used for the advective fluxes. This depends on the chosen discretization method.
Definition: flux/upwindscheme.hh:42
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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
Dune::FieldVector< Scalar, numComponents > molecularDiffusionFlux(const int phaseIdx) const
Returns the diffusive fluxes computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:129
Scalar heatConductionFlux(const int phaseIdx) const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:173
GetPropType< TypeTag, Properties::AdvectionType > AdvectionType
Definition: porousmediumflow/fluxvariables.hh:73
static constexpr bool enableAdvection
Definition: porousmediumflow/fluxvariables.hh:77
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition: porousmediumflow/fluxvariables.hh:74
Scalar advectiveFlux(const int phaseIdx, const FunctionType &upwindTerm) const
Returns the advective flux computed by the respective law. Specialization for enabled advection.
Definition: porousmediumflow/fluxvariables.hh:94
static constexpr bool enableEnergyBalance
Definition: porousmediumflow/fluxvariables.hh:79
Scalar heatConductionFlux() const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:157
PorousMediumFluxVariables()
The constructor.
Definition: porousmediumflow/fluxvariables.hh:83
static constexpr bool enableThermalNonEquilibrium
Definition: porousmediumflow/fluxvariables.hh:80
static constexpr bool enableMolecularDiffusion
Definition: porousmediumflow/fluxvariables.hh:78
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: porousmediumflow/fluxvariables.hh:75
Declares all properties used in Dumux.
Base class for the upwind scheme.