13#ifndef DUMUX_DISCRETIZATION_CC_TPFA_DISPERSION_FLUX_HH 
   14#define DUMUX_DISCRETIZATION_CC_TPFA_DISPERSION_FLUX_HH 
   16#include <dune/common/fvector.hh> 
   17#include <dune/common/fmatrix.hh> 
   30template<
class TypeTag, 
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
 
   37template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
 
   45    using FVElementGeometry = 
typename GridGeometry::LocalView;
 
   46    using SubControlVolume = 
typename GridGeometry::SubControlVolume;
 
   47    using SubControlVolumeFace = 
typename GridGeometry::SubControlVolumeFace;
 
   51    using ElementFluxVariablesCache = 
typename GridFluxVariablesCache::LocalView;
 
   52    using FluxVarCache = 
typename GridFluxVariablesCache::FluxVariablesCache;
 
   57    using Element = 
typename GridView::template Codim<0>::Entity;
 
   59    using Indices = 
typename ModelTraits::Indices;
 
   61    static constexpr int dim = GridView::dimension;
 
   62    static constexpr int dimWorld = GridView::dimensionworld;
 
   64    static constexpr int numPhases = ModelTraits::numFluidPhases();
 
   65    static constexpr int numComponents = ModelTraits::numFluidComponents();
 
   67    using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
 
   68    using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
 
   69    using HeatFluxScalar = Scalar;
 
   77    { 
return referenceSystem; }
 
 
   86                                                           const Element& element,
 
   87                                                           const FVElementGeometry& fvGeometry,
 
   88                                                           const ElementVolumeVariables& elemVolVars,
 
   89                                                           const SubControlVolumeFace& scvf,
 
   91                                                           const ElementFluxVariablesCache& elemFluxVarsCache)
 
   93        if (scvf.numOutsideScvs() > 1 )
 
   94            DUNE_THROW(Dune::NotImplemented, 
"\n Dispersion using ccTPFA is only implemented for conforming grids.");
 
   95        if (!stationaryVelocityField)
 
   96            DUNE_THROW(Dune::NotImplemented, 
"\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
 
   98        ComponentFluxVector componentFlux(0.0);
 
   99        const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
 
  100        const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
 
  103        const auto rhoOutside = 
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
 
  104        const Scalar rho = 0.5*(rhoInside + rhoOutside);
 
  106        for (
int compIdx = 0; compIdx < numComponents; compIdx++)
 
  108            if constexpr (!FluidSystem::isTracerFluidSystem())
 
  109                if (compIdx == FluidSystem::getMainComponent(phaseIdx))
 
  112            const auto& dispersionTensor =
 
  113                ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
 
  114                                                                                         elemVolVars, elemFluxVarsCache,
 
  119            const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor);
 
  121            const auto xInside = 
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
 
  122            const auto xOutide = 
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
 
  124            componentFlux[compIdx] = (rho * tij * (xInside-xOutide));
 
  125            if constexpr (!FluidSystem::isTracerFluidSystem())
 
  126                if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
 
  127                    componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
 
  129        return componentFlux;
 
 
  137                                                const Element& element,
 
  138                                                const FVElementGeometry& fvGeometry,
 
  139                                                const ElementVolumeVariables& elemVolVars,
 
  140                                                const SubControlVolumeFace& scvf,
 
  142                                                const ElementFluxVariablesCache& elemFluxVarsCache)
 
  144        if (scvf.numOutsideScvs() > 1 )
 
  145            DUNE_THROW(Dune::NotImplemented, 
"\n Dispersion using ccTPFA is only implemented for conforming grids.");
 
  146        if (!stationaryVelocityField)
 
  147            DUNE_THROW(Dune::NotImplemented, 
"\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
 
  149        const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
 
  150        const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
 
  152        const auto& dispersionTensor =
 
  153            ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
 
  154                                                                         elemVolVars, elemFluxVarsCache,
 
  159        const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor);
 
  162        const auto tInside = insideVolVars.temperature();
 
  163        const auto tOutside = outsideVolVars.temperature();
 
  166        return tij * (tInside-tOutside);
 
 
  172    template<
class Tensor>
 
  173    static Scalar calculateTransmissibility_(
const FVElementGeometry& fvGeometry,
 
  174                                             const ElementVolumeVariables& elemVolVars,
 
  175                                             const SubControlVolumeFace& scvf,
 
  180        const auto insideScvIdx = scvf.insideScvIdx();
 
  181        const auto& insideScv = fvGeometry.scv(insideScvIdx);
 
  182        const auto& insideVolVars = elemVolVars[insideScvIdx];
 
  189            tij = Extrusion::area(fvGeometry, scvf)*ti;
 
  194            const auto outsideScvIdx = scvf.outsideScvIdx();
 
  195            const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
 
  196            const auto& outsideVolVars = elemVolVars[outsideScvIdx];
 
  199            if constexpr (dim == dimWorld)
 
  209                tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
 
 
static ComponentFluxVector compositionalDispersionFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the dispersive fluxes of all components within a fluid phase across the given sub-control vol...
Definition cctpfa/dispersionflux.hh:85
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition cctpfa/dispersionflux.hh:76
static HeatFluxScalar thermalDispersionFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the thermal dispersive flux across the given sub-control volume face.
Definition cctpfa/dispersionflux.hh:136
Definition box/dispersionflux.hh:32
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition tpfa/computetransmissibility.hh:36
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition referencesystemformulation.hh:54
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition referencesystemformulation.hh:43
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition referencesystemformulation.hh:33
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:166
Traits of a flux variables type.
Definition flux/traits.hh:32
static constexpr bool hasStationaryVelocityField()
Definition flux/traits.hh:33
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...