43    using GridView = 
typename GridGeometry::GridView;
 
   44    using FVElementGeometry = 
typename GridGeometry::LocalView;
 
   45    using SubControlVolume = 
typename GridGeometry::SubControlVolume;
 
   46    using SubControlVolumeFace = 
typename GridGeometry::SubControlVolumeFace;
 
   48    using Element = 
typename GridView::template Codim<0>::Entity;
 
   50    static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
 
   51    static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
 
   54    static constexpr bool isSolDependent = (diffusionEnabled && diffusionIsSolDependent) ||
 
   55                                           (heatConductionEnabled && heatConductionIsSolDependent);
 
   59    : problemPtr_(&problem) {}
 
 
   72    template<
class FluxVariablesCacheContainer, 
class FluxVariablesCache, 
class ElementVolumeVariables>
 
   73    void fill(FluxVariablesCacheContainer& fluxVarsCacheContainer,
 
   74              FluxVariablesCache& scvfFluxVarsCache,
 
   75              const Element& element,
 
   76              const FVElementGeometry& fvGeometry,
 
   77              const ElementVolumeVariables& elemVolVars,
 
   78              const SubControlVolumeFace& scvf,
 
   79              const bool forceUpdateAll = 
false)
 
   84            if constexpr (diffusionEnabled)
 
   85                fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
 
   86            if constexpr (heatConductionEnabled)
 
   87                fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
 
   91            if constexpr (diffusionEnabled && diffusionIsSolDependent)
 
   92                fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
 
   93            if constexpr (heatConductionEnabled && heatConductionIsSolDependent)
 
   94                fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
 
 
  100    const Problem& problem()
 const 
  101    { 
return *problemPtr_; }
 
  105    template<
class FluxVariablesCache, 
class ElementVolumeVariables>
 
  106    void fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
 
  107                        const Element& element,
 
  108                        const FVElementGeometry& fvGeometry,
 
  109                        const ElementVolumeVariables& elemVolVars,
 
  110                        const SubControlVolumeFace& scvf)
 
  112        using DiffusionType = 
typename ElementVolumeVariables::VolumeVariables::MolecularDiffusionType;
 
  113        using DiffusionFiller = 
typename DiffusionType::Cache::Filler;
 
  114        using FluidSystem = 
typename ElementVolumeVariables::VolumeVariables::FluidSystem;
 
  116        static constexpr int numPhases = ModelTraits::numFluidPhases();
 
  117        static constexpr int numComponents = ModelTraits::numFluidComponents();
 
  120        if constexpr (FluidSystem::isTracerFluidSystem())
 
  121            for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
 
  122                for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
 
  123                    DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
 
  125            for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
 
  126                for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
 
  127                    if (compIdx != FluidSystem::getMainComponent(phaseIdx))
 
  128                        DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
 
  132    template<
class FluxVariablesCache, 
class ElementVolumeVariables>
 
  133    void fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
 
  134                             const Element& element,
 
  135                             const FVElementGeometry& fvGeometry,
 
  136                             const ElementVolumeVariables& elemVolVars,
 
  137                             const SubControlVolumeFace& scvf)
 
  139        using HeatConductionType = 
typename ElementVolumeVariables::VolumeVariables::HeatConductionType;
 
  140        using HeatConductionFiller = 
typename HeatConductionType::Cache::Filler;
 
  143        HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
 
  146    const Problem* problemPtr_;
 
 
void fill(FluxVariablesCacheContainer &fluxVarsCacheContainer, FluxVariablesCache &scvfFluxVarsCache, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const bool forceUpdateAll=false)
function to fill the flux variables caches
Definition scalarfluxvariablescachefiller.hh:73
FreeFlowScalarFluxVariablesCacheFillerImplementation< Problem, ModelTraits, diffusionIsSolDependent, heatConductionIsSolDependent, typename ProblemTraits< Problem >::GridGeometry::DiscretizationMethod > FreeFlowScalarFluxVariablesCacheFiller
The flux variables cache filler class for free flow.
Definition scalarfluxvariablescachefiller.hh:36