version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_POROUSMEDIUMFLOW_FLUXVARIABLES_HH
14#define DUMUX_POROUSMEDIUMFLOW_FLUXVARIABLES_HH
15
16#include <bitset>
17#include <array>
18
22
23namespace Dumux {
24
34template<class TypeTag,
35 class UpScheme = UpwindScheme<GetPropType<TypeTag, Properties::GridGeometry>> >
37: public FluxVariablesBase<GetPropType<TypeTag, Properties::Problem>,
38 typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView,
39 typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView,
40 typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView>
41{
44
45 enum
46 {
47 numPhases = ModelTraits::numFluidPhases(),
48 numComponents = ModelTraits::numFluidComponents()
49 };
50
51public:
52 using UpwindScheme = UpScheme;
57
58 static constexpr bool enableAdvection = ModelTraits::enableAdvection();
59 static constexpr bool enableMolecularDiffusion = ModelTraits::enableMolecularDiffusion();
60 static constexpr bool enableCompositionalDispersion = ModelTraits::enableCompositionalDispersion();
61 static constexpr bool enableThermalDispersion = ModelTraits::enableThermalDispersion();
62 static constexpr bool enableEnergyBalance = ModelTraits::enableEnergyBalance();
63 static constexpr bool enableThermalNonEquilibrium = getPropValue<TypeTag, Properties::EnableThermalNonEquilibrium>();
64
67 {
68 advFluxIsCached_.reset();
69 advFluxBeforeUpwinding_.fill(0.0);
70 }
71
75 template<typename FunctionType>
76 Scalar advectiveFlux([[maybe_unused]] const int phaseIdx, [[maybe_unused]] const FunctionType& upwindTerm) const
77 {
78 if constexpr (enableAdvection)
79 {
80 if (!advFluxIsCached_[phaseIdx])
81 {
82
83 advFluxBeforeUpwinding_[phaseIdx] = AdvectionType::flux(this->problem(),
84 this->element(),
85 this->fvGeometry(),
86 this->elemVolVars(),
87 this->scvFace(),
88 phaseIdx,
89 this->elemFluxVarsCache());
90 advFluxIsCached_.set(phaseIdx, true);
91 }
92
94 return UpwindScheme::apply(*this, upwindTerm, advFluxBeforeUpwinding_[phaseIdx], phaseIdx);
95 }
96 else
97 return 0.0;
98 }
99
103 Dune::FieldVector<Scalar, numComponents> molecularDiffusionFlux([[maybe_unused]] const int phaseIdx) const
104 {
105 if constexpr (enableMolecularDiffusion)
106 return MolecularDiffusionType::flux(this->problem(),
107 this->element(),
108 this->fvGeometry(),
109 this->elemVolVars(),
110 this->scvFace(),
111 phaseIdx,
112 this->elemFluxVarsCache());
113 else
114 return Dune::FieldVector<Scalar, numComponents>(0.0);
115 }
116
120 Dune::FieldVector<Scalar, numComponents> compositionalDispersionFlux([[maybe_unused]] const int phaseIdx) const
121 {
122 if constexpr (enableCompositionalDispersion)
123 {
124 return DispersionFluxType::compositionalDispersionFlux(this->problem(),
125 this->element(),
126 this->fvGeometry(),
127 this->elemVolVars(),
128 this->scvFace(),
129 phaseIdx,
130 this->elemFluxVarsCache());
131 }
132 else
133 return Dune::FieldVector<Scalar, numComponents>(0.0);
134 }
135
139 Dune::FieldVector<Scalar, 1> thermalDispersionFlux([[maybe_unused]] const int phaseIdx = 0) const
140 {
141 if constexpr (enableThermalDispersion)
142 {
143 return DispersionFluxType::thermalDispersionFlux(this->problem(),
144 this->element(),
145 this->fvGeometry(),
146 this->elemVolVars(),
147 this->scvFace(),
148 phaseIdx,
149 this->elemFluxVarsCache());
150 }
151 else
152 return Dune::FieldVector<Scalar, 1>(0.0);
153 }
154
159 Scalar heatConductionFlux() const
160 {
161 static_assert(!enableThermalNonEquilibrium, "This only works for thermal equilibrium");
162 if constexpr (enableEnergyBalance)
163 return HeatConductionType::flux(this->problem(),
164 this->element(),
165 this->fvGeometry(),
166 this->elemVolVars(),
167 this->scvFace(),
168 this->elemFluxVarsCache());
169 else
170 return 0.0;
171 }
172
177 Scalar heatConductionFlux([[maybe_unused]] const int phaseIdx) const
178 {
179 if constexpr (enableEnergyBalance)
180 return HeatConductionType::flux(this->problem(),
181 this->element(),
182 this->fvGeometry(),
183 this->elemVolVars(),
184 this->scvFace(),
185 phaseIdx,
186 this->elemFluxVarsCache());
187 else
188 return 0.0;
189 }
190
191private:
193 mutable std::bitset<numPhases> advFluxIsCached_;
194 mutable std::array<Scalar, numPhases> advFluxBeforeUpwinding_;
195};
196
197} // end namespace Dumux
198
199#endif
Base class for the flux variables living on a sub control volume face.
Definition: fluxvariablesbase.hh:33
The porous medium flux variables class that computes advective / convective, molecular diffusive and ...
Definition: porousmediumflow/fluxvariables.hh:41
Scalar advectiveFlux(const int phaseIdx, const FunctionType &upwindTerm) const
Returns the advective flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:76
Dune::FieldVector< Scalar, numComponents > compositionalDispersionFlux(const int phaseIdx) const
Returns the compositional dispersion flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:120
Dune::FieldVector< Scalar, numComponents > molecularDiffusionFlux(const int phaseIdx) const
Returns the diffusive fluxes computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:103
static constexpr bool enableAdvection
Definition: porousmediumflow/fluxvariables.hh:58
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: porousmediumflow/fluxvariables.hh:56
static constexpr bool enableThermalNonEquilibrium
Definition: porousmediumflow/fluxvariables.hh:63
static constexpr bool enableEnergyBalance
Definition: porousmediumflow/fluxvariables.hh:62
PorousMediumFluxVariables()
The constructor.
Definition: porousmediumflow/fluxvariables.hh:66
GetPropType< TypeTag, Properties::MolecularDiffusionType > MolecularDiffusionType
Definition: porousmediumflow/fluxvariables.hh:55
GetPropType< TypeTag, Properties::DispersionFluxType > DispersionFluxType
Definition: porousmediumflow/fluxvariables.hh:54
Dune::FieldVector< Scalar, 1 > thermalDispersionFlux(const int phaseIdx=0) const
Returns the thermal dispersion flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:139
static constexpr bool enableThermalDispersion
Definition: porousmediumflow/fluxvariables.hh:61
UpScheme UpwindScheme
Definition: porousmediumflow/fluxvariables.hh:52
Scalar heatConductionFlux(const int phaseIdx) const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:177
Scalar heatConductionFlux() const
Returns the conductive flux computed by the respective law.
Definition: porousmediumflow/fluxvariables.hh:159
GetPropType< TypeTag, Properties::AdvectionType > AdvectionType
Definition: porousmediumflow/fluxvariables.hh:53
static constexpr bool enableMolecularDiffusion
Definition: porousmediumflow/fluxvariables.hh:59
static constexpr bool enableCompositionalDispersion
Definition: porousmediumflow/fluxvariables.hh:60
Defines all properties used in Dumux.
Base class for the upwind scheme.
Base class for the flux variables living on a sub control volume face.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17