12#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH 
   13#define DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH 
   15#include <dune/common/fvector.hh> 
   30template<
class TypeTag, 
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
 
   37template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
 
   40    using Implementation = FicksLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>;
 
   44    using FVElementGeometry = 
typename GridGeometry::LocalView;
 
   45    using SubControlVolumeFace = 
typename GridGeometry::SubControlVolumeFace;
 
   47    using GridView = 
typename GridGeometry::GridView;
 
   48    using ElementVolumeVariables = 
typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
 
   49    using VolumeVariables = 
typename ElementVolumeVariables::VolumeVariables;
 
   50    using Element = 
typename GridView::template Codim<0>::Entity;
 
   52    using ElementFluxVariablesCache = 
typename GridFluxVariablesCache::LocalView;
 
   53    using FluxVariablesCache = 
typename GridFluxVariablesCache::FluxVariablesCache;
 
   57    static const int dim = GridView::dimension;
 
   58    static const int dimWorld = GridView::dimensionworld;
 
   59    static const int numPhases = ModelTraits::numFluidPhases();
 
   60    static const int numComponents = ModelTraits::numFluidComponents();
 
   62    using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
 
   65    class TpfaFicksLawCacheFiller
 
   70        template<
class FluxVariablesCacheFiller>
 
   71        static void fill(FluxVariablesCache& scvfFluxVarsCache,
 
   72                         unsigned int phaseIdx, 
unsigned int compIdx,
 
   73                         const Problem& problem,
 
   74                         const Element& element,
 
   75                         const FVElementGeometry& fvGeometry,
 
   76                         const ElementVolumeVariables& elemVolVars,
 
   77                         const SubControlVolumeFace& scvf,
 
   78                         const FluxVariablesCacheFiller& fluxVarsCacheFiller)
 
   80            scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
 
   85    class TpfaFicksLawCache
 
   88        using Filler = TpfaFicksLawCacheFiller;
 
   90        void updateDiffusion(
const Problem& problem,
 
   91                             const Element& element,
 
   92                             const FVElementGeometry& fvGeometry,
 
   93                             const ElementVolumeVariables& elemVolVars,
 
   94                             const SubControlVolumeFace &scvf,
 
   95                             const unsigned int phaseIdx,
 
   96                             const unsigned int compIdx)
 
   98            tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
 
  101        const Scalar& diffusionTij(
unsigned int phaseIdx, 
unsigned int compIdx)
 const 
  102        { 
return tij_[phaseIdx][compIdx]; }
 
  105        std::array< std::array<Scalar, numComponents>, numPhases> tij_;
 
  109    using DiscretizationMethod = DiscretizationMethods::CCTpfa;
 
  111    static constexpr DiscretizationMethod discMethod{};
 
  115    { 
return referenceSystem; }
 
  118    using Cache = TpfaFicksLawCache;
 
  120    using DiffusionCoefficientsContainer = FickianDiffusionCoefficients<Scalar, numPhases, numComponents>;
 
  132    static ComponentFluxVector flux(
const Problem& problem,
 
  133                                    const Element& element,
 
  134                                    const FVElementGeometry& fvGeometry,
 
  135                                    const VolumeVariables& insideVolVars,
 
  136                                    const VolumeVariables& outsideVolVars,
 
  137                                    const SubControlVolumeFace& scvf,
 
  139                                    const ElementFluxVariablesCache& elemFluxVarsCache)
 
  141        if constexpr (isMixedDimensional_)
 
  142            if (scvf.numOutsideScv() != 1)
 
  143                DUNE_THROW(Dune::Exception, 
"This flux overload requires scvf.numOutsideScv() == 1");
 
  146        const auto getOutsideX = [&](
const Scalar xInside, 
const Scalar tij, 
const int phaseIdx, 
const int compIdx)
 
  152        const auto getRho = [&](
const int phaseIdx, 
const Scalar rhoInside, 
const Scalar rhoOutside)
 
  154            return 0.5*(rhoInside + rhoOutside);
 
  157        return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
 
  167    static ComponentFluxVector flux(
const Problem& problem,
 
  168                                    const Element& element,
 
  169                                    const FVElementGeometry& fvGeometry,
 
  170                                    const ElementVolumeVariables& elemVolVars,
 
  171                                    const SubControlVolumeFace& scvf,
 
  173                                    const ElementFluxVariablesCache& elemFluxVarsCache)
 
  176        const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
 
  177        const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
 
  180        const auto getOutsideX = [&](
const Scalar xInside, 
const Scalar tij, 
const int phaseIdx, 
const int compIdx)
 
  182            const Scalar massOrMoleFractionOutside = 
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
 
  183            if constexpr (isMixedDimensional_)
 
  185                return scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
 
  186                                                : branchingFacetX_(problem, element, fvGeometry, elemVolVars,
 
  187                                                                   elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
 
  190                return massOrMoleFractionOutside;
 
  194        const auto getRho = [&](
const int phaseIdx, 
const Scalar rhoInside, 
const Scalar rhoOutside)
 
  196            if constexpr (isMixedDimensional_)
 
  198                return scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
 
  199                                                  : branchingFacetDensity_(elemVolVars, scvf, phaseIdx, rhoInside);
 
  202                return 0.5*(rhoInside + rhoOutside);
 
  205        return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
 
  209    static Scalar calculateTransmissibility(
const Problem& problem,
 
  210                                            const Element& element,
 
  211                                            const FVElementGeometry& fvGeometry,
 
  212                                            const ElementVolumeVariables& elemVolVars,
 
  213                                            const SubControlVolumeFace& scvf,
 
  214                                            const int phaseIdx, 
const int compIdx)
 
  218        const auto insideScvIdx = scvf.insideScvIdx();
 
  219        const auto& insideScv = fvGeometry.scv(insideScvIdx);
 
  220        const auto& insideVolVars = elemVolVars[insideScvIdx];
 
  221        const auto getDiffCoeff = [&](
const auto& vv)
 
  223            using FluidSystem = 
typename ElementVolumeVariables::VolumeVariables::FluidSystem;
 
  224            if constexpr (FluidSystem::isTracerFluidSystem())
 
  225                return vv.effectiveDiffusionCoefficient(0, 0, compIdx);
 
  227                return vv.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
 
  230        const auto insideDiffCoeff = getDiffCoeff(insideVolVars);
 
  236        if (scvf.boundary() || scvf.numOutsideScvs() > 1)
 
  237            tij = Extrusion::area(fvGeometry, scvf)*ti;
 
  242            const auto outsideScvIdx = scvf.outsideScvIdx();
 
  243            const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
 
  244            const auto& outsideVolVars = elemVolVars[outsideScvIdx];
 
  245            const auto outsideDiffCoeff = getDiffCoeff(outsideVolVars);
 
  248            if constexpr (dim == dimWorld)
 
  255                                                 outsideVolVars.extrusionFactor());
 
  261                tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
 
  268    template<
class Outs
ideFractionHelper, 
class DensityHelper>
 
  269    static ComponentFluxVector flux_(
const Element& element,
 
  270                                     const FVElementGeometry& fvGeometry,
 
  271                                     const VolumeVariables& insideVolVars,
 
  272                                     const VolumeVariables& outsideVolVars,
 
  273                                     const ElementFluxVariablesCache& elemFluxVarsCache,
 
  274                                     const SubControlVolumeFace& scvf,
 
  276                                     const OutsideFractionHelper& getOutsideX,
 
  277                                     const DensityHelper& getRho)
 
  279        ComponentFluxVector componentFlux(0.0);
 
  280        for (
int compIdx = 0; compIdx < numComponents; compIdx++)
 
  282            using FluidSystem = 
typename VolumeVariables::FluidSystem;
 
  283            if constexpr (!FluidSystem::isTracerFluidSystem())
 
  284                if (compIdx == FluidSystem::getMainComponent(phaseIdx))
 
  288            const Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
 
  291            const Scalar xInside = 
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
 
  292            const Scalar xOutside = getOutsideX(xInside, tij, phaseIdx, compIdx);
 
  294            const Scalar rhoInside = 
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
 
  295            const Scalar rhoOutside = 
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
 
  297            const Scalar rho = getRho(phaseIdx, rhoInside, rhoOutside);
 
  299            componentFlux[compIdx] = rho*tij*(xInside - xOutside);
 
  300            if constexpr (!FluidSystem::isTracerFluidSystem())
 
  301                if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
 
  302                    componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
 
  305        return componentFlux;
 
  309    static Scalar branchingFacetX_(
const Problem& problem,
 
  310                                   const Element& element,
 
  311                                   const FVElementGeometry& fvGeometry,
 
  312                                   const ElementVolumeVariables& elemVolVars,
 
  313                                   const ElementFluxVariablesCache& elemFluxVarsCache,
 
  314                                   const SubControlVolumeFace& scvf,
 
  315                                   const Scalar insideX, 
const Scalar insideTi,
 
  316                                   const int phaseIdx, 
const int compIdx)
 
  318        Scalar sumTi(insideTi);
 
  319        Scalar sumXTi(insideTi*insideX);
 
  321        for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
 
  323            const auto outsideScvIdx = scvf.outsideScvIdx(i);
 
  324            const auto& outsideVolVars = elemVolVars[outsideScvIdx];
 
  325            const Scalar massOrMoleFractionOutside = 
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
 
  326            const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
 
  328            const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
 
  330            sumXTi += outsideTi*massOrMoleFractionOutside;
 
  333        return sumTi > 0 ? sumXTi/sumTi : 0;
 
  337    static Scalar branchingFacetDensity_(
const ElementVolumeVariables& elemVolVars,
 
  338                                         const SubControlVolumeFace& scvf,
 
  340                                         const Scalar insideRho)
 
  342        Scalar rho(insideRho);
 
  343        for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
 
  345            const auto outsideScvIdx = scvf.outsideScvIdx(i);
 
  346            const auto& outsideVolVars = elemVolVars[outsideScvIdx];
 
  347            const Scalar rhoOutside = 
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
 
  350        return rho/(scvf.numOutsideScvs()+1);
 
  353    static constexpr bool isMixedDimensional_ = 
static_cast<int>(GridView::dimension) < 
static_cast<int>(GridView::dimensionworld);
 
forward declaration of the method-specific implementation
Definition flux/box/fickslaw.hh:32
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
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
The available discretization methods in Dumux.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...