3.2-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
92 template<typename FunctionType>
93 Scalar advectiveFlux([[maybe_unused]] const int phaseIdx, [[maybe_unused]] const FunctionType& upwindTerm) const
94 {
95 if constexpr (enableAdvection)
96 {
97 if (!advFluxIsCached_[phaseIdx])
98 {
99
100 advFluxBeforeUpwinding_[phaseIdx] = AdvectionType::flux(this->problem(),
101 this->element(),
102 this->fvGeometry(),
103 this->elemVolVars(),
104 this->scvFace(),
105 phaseIdx,
106 this->elemFluxVarsCache());
107 advFluxIsCached_.set(phaseIdx, true);
108 }
109
111 return UpwindScheme::apply(*this, upwindTerm, advFluxBeforeUpwinding_[phaseIdx], phaseIdx);
112 }
113 else
114 return 0.0;
115 }
116
120 Dune::FieldVector<Scalar, numComponents> molecularDiffusionFlux([[maybe_unused]] const int phaseIdx) const
121 {
122 if constexpr (enableMolecularDiffusion)
123 return MolecularDiffusionType::flux(this->problem(),
124 this->element(),
125 this->fvGeometry(),
126 this->elemVolVars(),
127 this->scvFace(),
128 phaseIdx,
129 this->elemFluxVarsCache());
130 else
131 return {0.0};
132 }
133
138 Scalar heatConductionFlux() const
139 {
140 static_assert(!enableThermalNonEquilibrium, "This only works for thermal equilibrium");
141 if constexpr (enableEnergyBalance)
142 return HeatConductionType::flux(this->problem(),
143 this->element(),
144 this->fvGeometry(),
145 this->elemVolVars(),
146 this->scvFace(),
147 this->elemFluxVarsCache());
148 else
149 return 0.0;
150 }
151
156 Scalar heatConductionFlux([[maybe_unused]] const int phaseIdx) const
157 {
158 if constexpr (enableEnergyBalance)
159 return HeatConductionType::flux(this->problem(),
160 this->element(),
161 this->fvGeometry(),
162 this->elemVolVars(),
163 this->scvFace(),
164 phaseIdx,
165 this->elemFluxVarsCache());
166 else
167 return 0.0;
168 }
169
170private:
172 mutable std::bitset<numPhases> advFluxIsCached_;
173 mutable std::array<Scalar, numPhases> advFluxBeforeUpwinding_;
174};
175
176} // end namespace Dumux
177
178#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:156
GetPropType< TypeTag, Properties::AdvectionType > AdvectionType
Definition: porousmediumflow/fluxvariables.hh:73
static constexpr bool enableAdvection
Definition: porousmediumflow/fluxvariables.hh:77
Scalar advectiveFlux(const int phaseIdx, const FunctionType &upwindTerm) const
Returns the advective flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:93
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition: porousmediumflow/fluxvariables.hh:74
Scalar heatConductionFlux() const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:138
static constexpr bool enableEnergyBalance
Definition: porousmediumflow/fluxvariables.hh:79
Dune::FieldVector< Scalar, numComponents > molecularDiffusionFlux(const int phaseIdx) const
Returns the diffusive fluxes computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:120
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.