3.3.0
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{
56
57 enum
58 {
59 numPhases = ModelTraits::numFluidPhases(),
60 numComponents = ModelTraits::numFluidComponents()
61 };
62
63public:
67
68 static constexpr bool enableAdvection = ModelTraits::enableAdvection();
69 static constexpr bool enableMolecularDiffusion = ModelTraits::enableMolecularDiffusion();
70 static constexpr bool enableEnergyBalance = ModelTraits::enableEnergyBalance();
71 static constexpr bool enableThermalNonEquilibrium = getPropValue<TypeTag, Properties::EnableThermalNonEquilibrium>();
72
75 {
76 advFluxIsCached_.reset();
77 advFluxBeforeUpwinding_.fill(0.0);
78 }
79
83 template<typename FunctionType>
84 Scalar advectiveFlux([[maybe_unused]] const int phaseIdx, [[maybe_unused]] const FunctionType& upwindTerm) const
85 {
86 if constexpr (enableAdvection)
87 {
88 if (!advFluxIsCached_[phaseIdx])
89 {
90
91 advFluxBeforeUpwinding_[phaseIdx] = AdvectionType::flux(this->problem(),
92 this->element(),
93 this->fvGeometry(),
94 this->elemVolVars(),
95 this->scvFace(),
96 phaseIdx,
97 this->elemFluxVarsCache());
98 advFluxIsCached_.set(phaseIdx, true);
99 }
100
102 return UpwindScheme::apply(*this, upwindTerm, advFluxBeforeUpwinding_[phaseIdx], phaseIdx);
103 }
104 else
105 return 0.0;
106 }
107
111 Dune::FieldVector<Scalar, numComponents> molecularDiffusionFlux([[maybe_unused]] const int phaseIdx) const
112 {
113 if constexpr (enableMolecularDiffusion)
114 return MolecularDiffusionType::flux(this->problem(),
115 this->element(),
116 this->fvGeometry(),
117 this->elemVolVars(),
118 this->scvFace(),
119 phaseIdx,
120 this->elemFluxVarsCache());
121 else
122 return Dune::FieldVector<Scalar, numComponents>(0.0);
123 }
124
129 Scalar heatConductionFlux() const
130 {
131 static_assert(!enableThermalNonEquilibrium, "This only works for thermal equilibrium");
132 if constexpr (enableEnergyBalance)
133 return HeatConductionType::flux(this->problem(),
134 this->element(),
135 this->fvGeometry(),
136 this->elemVolVars(),
137 this->scvFace(),
138 this->elemFluxVarsCache());
139 else
140 return 0.0;
141 }
142
147 Scalar heatConductionFlux([[maybe_unused]] const int phaseIdx) const
148 {
149 if constexpr (enableEnergyBalance)
150 return HeatConductionType::flux(this->problem(),
151 this->element(),
152 this->fvGeometry(),
153 this->elemVolVars(),
154 this->scvFace(),
155 phaseIdx,
156 this->elemFluxVarsCache());
157 else
158 return 0.0;
159 }
160
161private:
163 mutable std::bitset<numPhases> advFluxIsCached_;
164 mutable std::array<Scalar, numPhases> advFluxBeforeUpwinding_;
165};
166
167} // end namespace Dumux
168
169#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
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
Scalar heatConductionFlux(const int phaseIdx) const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:147
GetPropType< TypeTag, Properties::AdvectionType > AdvectionType
Definition: porousmediumflow/fluxvariables.hh:64
static constexpr bool enableAdvection
Definition: porousmediumflow/fluxvariables.hh:68
Scalar advectiveFlux(const int phaseIdx, const FunctionType &upwindTerm) const
Returns the advective flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:84
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition: porousmediumflow/fluxvariables.hh:65
Scalar heatConductionFlux() const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:129
static constexpr bool enableEnergyBalance
Definition: porousmediumflow/fluxvariables.hh:70
Dune::FieldVector< Scalar, numComponents > molecularDiffusionFlux(const int phaseIdx) const
Returns the diffusive fluxes computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:111
PorousMediumFluxVariables()
The constructor.
Definition: porousmediumflow/fluxvariables.hh:74
static constexpr bool enableThermalNonEquilibrium
Definition: porousmediumflow/fluxvariables.hh:71
static constexpr bool enableMolecularDiffusion
Definition: porousmediumflow/fluxvariables.hh:69
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: porousmediumflow/fluxvariables.hh:66
Declares all properties used in Dumux.
Base class for the upwind scheme.