3.6-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 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;
69
70 static constexpr bool enableAdvection = ModelTraits::enableAdvection();
71 static constexpr bool enableMolecularDiffusion = ModelTraits::enableMolecularDiffusion();
72 static constexpr bool enableCompositionalDispersion = ModelTraits::enableCompositionalDispersion();
73 static constexpr bool enableThermalDispersion = ModelTraits::enableThermalDispersion();
74 static constexpr bool enableEnergyBalance = ModelTraits::enableEnergyBalance();
75 static constexpr bool enableThermalNonEquilibrium = getPropValue<TypeTag, Properties::EnableThermalNonEquilibrium>();
76
79 {
80 advFluxIsCached_.reset();
81 advFluxBeforeUpwinding_.fill(0.0);
82 }
83
87 template<typename FunctionType>
88 Scalar advectiveFlux([[maybe_unused]] const int phaseIdx, [[maybe_unused]] const FunctionType& upwindTerm) const
89 {
90 if constexpr (enableAdvection)
91 {
92 if (!advFluxIsCached_[phaseIdx])
93 {
94
95 advFluxBeforeUpwinding_[phaseIdx] = AdvectionType::flux(this->problem(),
96 this->element(),
97 this->fvGeometry(),
98 this->elemVolVars(),
99 this->scvFace(),
100 phaseIdx,
101 this->elemFluxVarsCache());
102 advFluxIsCached_.set(phaseIdx, true);
103 }
104
106 return UpwindScheme::apply(*this, upwindTerm, advFluxBeforeUpwinding_[phaseIdx], phaseIdx);
107 }
108 else
109 return 0.0;
110 }
111
115 Dune::FieldVector<Scalar, numComponents> molecularDiffusionFlux([[maybe_unused]] const int phaseIdx) const
116 {
117 if constexpr (enableMolecularDiffusion)
118 return MolecularDiffusionType::flux(this->problem(),
119 this->element(),
120 this->fvGeometry(),
121 this->elemVolVars(),
122 this->scvFace(),
123 phaseIdx,
124 this->elemFluxVarsCache());
125 else
126 return Dune::FieldVector<Scalar, numComponents>(0.0);
127 }
128
132 Dune::FieldVector<Scalar, numComponents> compositionalDispersionFlux([[maybe_unused]] const int phaseIdx) const
133 {
134 if constexpr (enableCompositionalDispersion)
135 {
136 return DispersionFluxType::compositionalDispersionFlux(this->problem(),
137 this->element(),
138 this->fvGeometry(),
139 this->elemVolVars(),
140 this->scvFace(),
141 phaseIdx,
142 this->elemFluxVarsCache());
143 }
144 else
145 return Dune::FieldVector<Scalar, numComponents>(0.0);
146 }
147
151 Dune::FieldVector<Scalar, 1> thermalDispersionFlux([[maybe_unused]] const int phaseIdx = 0) const
152 {
153 if constexpr (enableThermalDispersion)
154 {
155 return DispersionFluxType::thermalDispersionFlux(this->problem(),
156 this->element(),
157 this->fvGeometry(),
158 this->elemVolVars(),
159 this->scvFace(),
160 phaseIdx,
161 this->elemFluxVarsCache());
162 }
163 else
164 return Dune::FieldVector<Scalar, 1>(0.0);
165 }
166
171 Scalar heatConductionFlux() const
172 {
173 static_assert(!enableThermalNonEquilibrium, "This only works for thermal equilibrium");
174 if constexpr (enableEnergyBalance)
175 return HeatConductionType::flux(this->problem(),
176 this->element(),
177 this->fvGeometry(),
178 this->elemVolVars(),
179 this->scvFace(),
180 this->elemFluxVarsCache());
181 else
182 return 0.0;
183 }
184
189 Scalar heatConductionFlux([[maybe_unused]] const int phaseIdx) const
190 {
191 if constexpr (enableEnergyBalance)
192 return HeatConductionType::flux(this->problem(),
193 this->element(),
194 this->fvGeometry(),
195 this->elemVolVars(),
196 this->scvFace(),
197 phaseIdx,
198 this->elemFluxVarsCache());
199 else
200 return 0.0;
201 }
202
203private:
205 mutable std::bitset<numPhases> advFluxIsCached_;
206 mutable std::array<Scalar, numPhases> advFluxBeforeUpwinding_;
207};
208
209} // end namespace Dumux
210
211#endif
Base class for the flux variables living on a sub control volume face.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
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:88
Dune::FieldVector< Scalar, numComponents > compositionalDispersionFlux(const int phaseIdx) const
Returns the compositional dispersion flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:132
Dune::FieldVector< Scalar, numComponents > molecularDiffusionFlux(const int phaseIdx) const
Returns the diffusive fluxes computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:115
static constexpr bool enableAdvection
Definition: porousmediumflow/fluxvariables.hh:70
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: porousmediumflow/fluxvariables.hh:68
static constexpr bool enableThermalNonEquilibrium
Definition: porousmediumflow/fluxvariables.hh:75
static constexpr bool enableEnergyBalance
Definition: porousmediumflow/fluxvariables.hh:74
PorousMediumFluxVariables()
The constructor.
Definition: porousmediumflow/fluxvariables.hh:78
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition: porousmediumflow/fluxvariables.hh:67
GetPropType< TypeTag, Properties::DispersionFluxType > DispersionFluxType
Definition: porousmediumflow/fluxvariables.hh:66
Dune::FieldVector< Scalar, 1 > thermalDispersionFlux(const int phaseIdx=0) const
Returns the thermal dispersion flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:151
static constexpr bool enableThermalDispersion
Definition: porousmediumflow/fluxvariables.hh:73
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:189
Scalar heatConductionFlux() const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:171
GetPropType< TypeTag, Properties::AdvectionType > AdvectionType
Definition: porousmediumflow/fluxvariables.hh:65
static constexpr bool enableMolecularDiffusion
Definition: porousmediumflow/fluxvariables.hh:71
static constexpr bool enableCompositionalDispersion
Definition: porousmediumflow/fluxvariables.hh:72
Declares all properties used in Dumux.
Base class for the upwind scheme.