12#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FICKS_LAW_HH 
   13#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FICKS_LAW_HH 
   18#include <dune/common/float_cmp.hh> 
   19#include <dune/common/fvector.hh> 
   35template<
class TypeTag,
 
   46template<
class TypeTag, ReferenceSystemFormulation referenceSystem =  ReferenceSystemFormulation::massAveraged>
 
   56template<
class TypeTag, ReferenceSystemFormulation referenceSystem>
 
   64    using FVElementGeometry = 
typename GridGeometry::LocalView;
 
   65    using SubControlVolume = 
typename GridGeometry::SubControlVolume;
 
   66    using SubControlVolumeFace = 
typename GridGeometry::SubControlVolumeFace;
 
   69    using GridView = 
typename GridGeometry::GridView;
 
   70    using Element = 
typename GridView::template Codim<0>::Entity;
 
   76    static const int numPhases = ModelTraits::numFluidPhases();
 
   77    static const int numComponents = ModelTraits::numFluidComponents();
 
   80    using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
 
   85    class FacetCouplingFicksLawCache
 
   89        using Filler = 
typename ParentType::Cache::Filler;
 
   94        static constexpr int insideTijIdx = 0;
 
   95        static constexpr int outsideTijIdx = 1;
 
   96        static constexpr int facetTijIdx = 2;
 
   99        using DiffusionTransmissibilityContainer = std::array<Scalar, 3>;
 
  102        template< 
class Problem, 
class ElementVolumeVariables >
 
  103        void updateDiffusion(
const Problem& problem,
 
  104                             const Element& element,
 
  105                             const FVElementGeometry& fvGeometry,
 
  106                             const ElementVolumeVariables& elemVolVars,
 
  107                             const SubControlVolumeFace &scvf,
 
  108                             unsigned int phaseIdx,
 
  109                             unsigned int compIdx)
 
  111            tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
 
  118        Scalar diffusionTij(
unsigned int phaseIdx, 
unsigned int compIdx)
 const 
  119        { 
return tij_[phaseIdx][compIdx][insideTijIdx]; }
 
  122        Scalar diffusionTijInside(
unsigned int phaseIdx, 
unsigned int compIdx)
 const 
  123        { 
return tij_[phaseIdx][compIdx][insideTijIdx]; }
 
  126        Scalar diffusionTijOutside(
unsigned int phaseIdx, 
unsigned int compIdx)
 const 
  127        { 
return tij_[phaseIdx][compIdx][outsideTijIdx]; }
 
  130        Scalar diffusionTijFacet(
unsigned int phaseIdx, 
unsigned int compIdx)
 const 
  131        { 
return tij_[phaseIdx][compIdx][facetTijIdx]; }
 
  134        std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
 
  139    using Cache = FacetCouplingFicksLawCache;
 
  143    static constexpr DiscretizationMethod discMethod{};
 
  147    { 
return referenceSystem; }
 
  150    template< 
class Problem, 
class ElementVolumeVariables, 
class ElementFluxVarsCache >
 
  151    static ComponentFluxVector flux(
const Problem& problem,
 
  152                                    const Element& element,
 
  153                                    const FVElementGeometry& fvGeometry,
 
  154                                    const ElementVolumeVariables& elemVolVars,
 
  155                                    const SubControlVolumeFace& scvf,
 
  157                                    const ElementFluxVarsCache& elemFluxVarsCache)
 
  159        if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
 
  160            return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
 
  162        ComponentFluxVector componentFlux(0.0);
 
  163        for (
int compIdx = 0; compIdx < numComponents; compIdx++)
 
  165            if(compIdx == FluidSystem::getMainComponent(phaseIdx))
 
  169            const auto& fluxVarsCache = elemFluxVarsCache[scvf];
 
  170            const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
 
  171            const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
 
  172            const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
 
  175            const Scalar xInside = 
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
 
  176            const Scalar xFacet = 
massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
 
  177            const Scalar xOutside = 
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
 
  179            const Scalar rhoInside = 
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
 
  181            const Scalar rho = 0.5*(rhoInside + rhoFacet);
 
  183            componentFlux[compIdx] = fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
 
  184                                     + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet;
 
  186            if (!scvf.boundary())
 
  187                componentFlux[compIdx] += fluxVarsCache.diffusionTijOutside(phaseIdx, compIdx)*xOutside;
 
  188            componentFlux[compIdx] *= rho;
 
  190            if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
 
  191                componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
 
  194        return componentFlux;
 
  199    template< 
class Problem, 
class ElementVolumeVariables >
 
  200    static typename Cache::DiffusionTransmissibilityContainer
 
  201    calculateTransmissibility(
const Problem& problem,
 
  202                              const Element& element,
 
  203                              const FVElementGeometry& fvGeometry,
 
  204                              const ElementVolumeVariables& elemVolVars,
 
  205                              const SubControlVolumeFace& scvf,
 
  206                              unsigned int phaseIdx, 
unsigned int compIdx)
 
  208        typename Cache::DiffusionTransmissibilityContainer tij;
 
  209        if (!problem.couplingManager().isCoupled(element, scvf))
 
  212            tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
 
  219        const auto insideScvIdx = scvf.insideScvIdx();
 
  220        const auto& insideScv = fvGeometry.scv(insideScvIdx);
 
  221        const auto& insideVolVars = elemVolVars[insideScvIdx];
 
  222        const auto wIn = Extrusion::area(fvGeometry, scvf)
 
  224                                                      insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
 
  225                                                      insideVolVars.extrusionFactor());
 
  228        const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
 
  231        if (iBcTypes.hasOnlyNeumann())
 
  233            const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
 
  234            const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
 
  235                                   /facetVolVars.extrusionFactor()
 
  236                                   *
vtmv(scvf.unitOuterNormal(),
 
  237                                         facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
 
  238                                         scvf.unitOuterNormal());
 
  247            if (!scvf.boundary())
 
  249                const auto outsideScvIdx = scvf.outsideScvIdx();
 
  250                const auto& outsideVolVars = elemVolVars[outsideScvIdx];
 
  251                const auto wOut = -1.0*Extrusion::area(fvGeometry, scvf)
 
  253                                                               outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
 
  254                                                               outsideVolVars.extrusionFactor());
 
  256                if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
 
  260                    const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
 
  261                    tij[Cache::insideTijIdx]  = factor * ( wOut * xi + wFacet );
 
  262                    tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
 
  263                    tij[Cache::facetTijIdx]   = factor * ( - wOut - wFacet );
 
  267                    tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
 
  268                    tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
 
  269                    tij[Cache::outsideTijIdx] = 0.0;
 
  274                tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
 
  275                tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
 
  276                tij[Cache::outsideTijIdx] = 0.0;
 
  279        else if (iBcTypes.hasOnlyDirichlet())
 
  281            tij[Cache::insideTijIdx] = wIn;
 
  282            tij[Cache::outsideTijIdx] = 0.0;
 
  283            tij[Cache::facetTijIdx] = -wIn;
 
  286            DUNE_THROW(Dune::NotImplemented, 
"Interior boundary types other than pure Dirichlet or Neumann");
 
  296template<
class TypeTag, ReferenceSystemFormulation referenceSystem>
 
  297class CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, true>
 
  298: 
public FicksLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>
 
  300    using Implementation = CCTpfaFacetCouplingFicksLawImpl<TypeTag, referenceSystem, true>;
 
  301    using ParentType = FicksLawImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>;
 
  303    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
 
  304    using FVElementGeometry = 
typename GridGeometry::LocalView;
 
  305    using SubControlVolume = 
typename GridGeometry::SubControlVolume;
 
  306    using SubControlVolumeFace = 
typename GridGeometry::SubControlVolumeFace;
 
  309    using GridView = 
typename GridGeometry::GridView;
 
  310    using Element = 
typename GridView::template Codim<0>::Entity;
 
  312    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
 
  313    using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
 
  314    using BalanceEqOpts = GetPropType<TypeTag, Properties::BalanceEqOpts>;
 
  316    static const int numPhases = ModelTraits::numFluidPhases();
 
  317    static const int numComponents = ModelTraits::numFluidComponents();
 
  319    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
 
  320    using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
 
  325    class FacetCouplingFicksLawCache
 
  329        using Filler = 
typename ParentType::Cache::Filler;
 
  334        static constexpr int insideTijIdx = 0;
 
  335        static constexpr int facetTijIdx = 1;
 
  338        using DiffusionTransmissibilityContainer = std::array<Scalar, 2>;
 
  341        template< 
class Problem, 
class ElementVolumeVariables >
 
  342        void updateDiffusion(
const Problem& problem,
 
  343                             const Element& element,
 
  344                             const FVElementGeometry& fvGeometry,
 
  345                             const ElementVolumeVariables& elemVolVars,
 
  346                             const SubControlVolumeFace &scvf,
 
  347                             unsigned int phaseIdx,
 
  348                             unsigned int compIdx)
 
  350            tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
 
  357        Scalar diffusionTij(
unsigned int phaseIdx, 
unsigned int compIdx)
 const 
  358        { 
return tij_[phaseIdx][compIdx][insideTijIdx]; }
 
  361        Scalar diffusionTijInside(
unsigned int phaseIdx, 
unsigned int compIdx)
 const 
  362        { 
return tij_[phaseIdx][compIdx][insideTijIdx]; }
 
  365        Scalar diffusionTijFacet(
unsigned int phaseIdx, 
unsigned int compIdx)
 const 
  366        { 
return tij_[phaseIdx][compIdx][facetTijIdx]; }
 
  369        std::array< std::array<DiffusionTransmissibilityContainer, numComponents>, numPhases> tij_;
 
  374    using Cache = FacetCouplingFicksLawCache;
 
  376    using DiscretizationMethod = DiscretizationMethods::CCTpfa;
 
  378    static constexpr DiscretizationMethod discMethod{};
 
  382    { 
return referenceSystem; }
 
  385    template< 
class Problem, 
class ElementVolumeVariables, 
class ElementFluxVarsCache >
 
  386    static ComponentFluxVector flux(
const Problem& problem,
 
  387                                    const Element& element,
 
  388                                    const FVElementGeometry& fvGeometry,
 
  389                                    const ElementVolumeVariables& elemVolVars,
 
  390                                    const SubControlVolumeFace& scvf,
 
  392                                    const ElementFluxVarsCache& elemFluxVarsCache)
 
  394        if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
 
  395            return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
 
  397        ComponentFluxVector componentFlux(0.0);
 
  398        for (
int compIdx = 0; compIdx < numComponents; compIdx++)
 
  400            if(compIdx == FluidSystem::getMainComponent(phaseIdx))
 
  404            const auto& fluxVarsCache = elemFluxVarsCache[scvf];
 
  405            const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
 
  406            const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
 
  409            const Scalar xInside = 
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
 
  410            const Scalar xFacet = 
massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx);
 
  412            const Scalar rhoInside = 
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
 
  414            const Scalar rho = 0.5*(rhoInside + rhoFacet);
 
  416            componentFlux[compIdx] = rho*(fluxVarsCache.diffusionTijInside(phaseIdx, compIdx)*xInside
 
  417                                          + fluxVarsCache.diffusionTijFacet(phaseIdx, compIdx)*xFacet);
 
  419            if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
 
  420                componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
 
  423        return componentFlux;
 
  428    template< 
class Problem, 
class ElementVolumeVariables >
 
  429    static typename Cache::DiffusionTransmissibilityContainer
 
  430    calculateTransmissibility(
const Problem& problem,
 
  431                              const Element& element,
 
  432                              const FVElementGeometry& fvGeometry,
 
  433                              const ElementVolumeVariables& elemVolVars,
 
  434                              const SubControlVolumeFace& scvf,
 
  435                              unsigned int phaseIdx, 
unsigned int compIdx)
 
  437        typename Cache::DiffusionTransmissibilityContainer tij;
 
  438        if (!problem.couplingManager().isCoupled(element, scvf))
 
  441            tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
 
  451        if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
 
  452            DUNE_THROW(Dune::InvalidStateException, 
"Xi != 1.0 cannot be used on surface grids");
 
  454        const auto insideScvIdx = scvf.insideScvIdx();
 
  455        const auto& insideScv = fvGeometry.scv(insideScvIdx);
 
  456        const auto& insideVolVars = elemVolVars[insideScvIdx];
 
  457        const auto wIn = Extrusion::area(fvGeometry, scvf)
 
  459                                                      insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
 
  460                                                      insideVolVars.extrusionFactor());
 
  463        const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
 
  466        if (iBcTypes.hasOnlyNeumann())
 
  471            const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
 
  472            const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
 
  473                                   /sqrt(facetVolVars.extrusionFactor())
 
  474                                   *
vtmv(scvf.unitOuterNormal(),
 
  475                                         facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx),
 
  476                                         scvf.unitOuterNormal());
 
  478            tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
 
  479            tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
 
  481        else if (iBcTypes.hasOnlyDirichlet())
 
  483            tij[Cache::insideTijIdx] = wIn;
 
  484            tij[Cache::facetTijIdx] = -wIn;
 
  487            DUNE_THROW(Dune::NotImplemented, 
"Interior boundary types other than pure Dirichlet or Neumann");
 
Forward declaration of the implementation.
Definition multidomain/facet/cellcentered/tpfa/fickslaw.hh:38
Traits extracting the public Extrusion type from T Defaults to NoExtrusion if no such type is found.
Definition extrusion.hh:155
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.
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
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
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition math.hh:880
CCTpfaFacetCouplingFicksLawImpl< TypeTag, referenceSystem,(int(GetPropType< TypeTag, Properties::GridGeometry >::GridView::dimension)< int(GetPropType< TypeTag, Properties::GridGeometry >::GridView::dimensionworld)) > CCTpfaFacetCouplingFicksLaw
Fick's law for cell-centered finite volume schemes with two-point flux approximation in the context o...
Definition multidomain/facet/cellcentered/tpfa/fickslaw.hh:47
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
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
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
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...