13#ifndef DUMUX_POROUSMEDIUMFLOW_FLUXVARIABLES_HH 
   14#define DUMUX_POROUSMEDIUMFLOW_FLUXVARIABLES_HH 
   34template<
class TypeTag,
 
   38                           typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView,
 
   39                           typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView,
 
   40                           typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView>
 
   47        numPhases = ModelTraits::numFluidPhases(),
 
   48        numComponents = ModelTraits::numFluidComponents()
 
   68        advFluxIsCached_.reset();
 
   69        advFluxBeforeUpwinding_.fill(0.0);
 
 
   75    template<
typename FunctionType>
 
   76    Scalar 
advectiveFlux([[maybe_unused]] 
const int phaseIdx, [[maybe_unused]] 
const FunctionType& upwindTerm)
 const 
   80            if (!advFluxIsCached_[phaseIdx])
 
   83                advFluxBeforeUpwinding_[phaseIdx] = AdvectionType::flux(this->
problem(),
 
   90                advFluxIsCached_.set(phaseIdx, 
true);
 
 
  106            return MolecularDiffusionType::flux(this->
problem(),
 
  114            return Dune::FieldVector<Scalar, numComponents>(0.0);
 
 
  124            return DispersionFluxType::compositionalDispersionFlux(this->
problem(),
 
  133            return Dune::FieldVector<Scalar, numComponents>(0.0);
 
 
  143            return DispersionFluxType::thermalDispersionFlux(this->
problem(),
 
  152            return Dune::FieldVector<Scalar, 1>(0.0);
 
 
  163            return HeatConductionType::flux(this->
problem(),
 
 
  180            return HeatConductionType::flux(this->
problem(),
 
 
  193    mutable std::bitset<numPhases> advFluxIsCached_;
 
  194    mutable std::array<Scalar, numPhases> advFluxBeforeUpwinding_;
 
 
static Scalar apply(const FluxVariables &fluxVars, const UpwindTermFunction &upwindTerm, Scalar flux, int phaseIdx)
Definition multidomain/facet/box/upwindscheme.hh:33
Base class for the flux variables living on a sub control volume face.
Definition fluxvariablesbase.hh:33
const GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView & elemVolVars() const
Definition fluxvariablesbase.hh:69
const SubControlVolumeFace & scvFace() const
Definition fluxvariablesbase.hh:63
const GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView & elemFluxVarsCache() const
Definition fluxvariablesbase.hh:72
const GetPropType< TypeTag, Properties::GridGeometry >::LocalView & fvGeometry() const
Definition fluxvariablesbase.hh:66
const Element & element() const
Definition fluxvariablesbase.hh:60
const Problem & problem() const
Definition fluxvariablesbase.hh:57
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
BoxFacetCouplingUpwindScheme< GetPropType< TypeTag, Properties::GridGeometry > > 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.
UpwindSchemeImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > UpwindScheme
The upwind scheme used for the advective fluxes. This depends on the chosen discretization method.
Definition flux/upwindscheme.hh:30
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:310
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296