24#ifndef DUMUX_FREEFLOW_SCALAR_FLUXVARIABLESCACHE_FILLER_HH
25#define DUMUX_FREEFLOW_SCALAR_FLUXVARIABLESCACHE_FILLER_HH
38template<
class Problem,
class ModelTraits,
bool diffusionIsSolDependent,
bool heatConductionIsSolDependent,
class DiscretizationMethod>
47template<
class Problem,
class ModelTraits,
bool diffusionIsSolDependent,
bool heatConductionIsSolDependent>
51template<
class Problem,
class ModelTraits,
bool diffusionIsSolDependent,
bool heatConductionIsSolDependent>
55 using GridView =
typename GridGeometry::GridView;
56 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolume =
typename GridGeometry::SubControlVolume;
58 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
60 using Element =
typename GridView::template Codim<0>::Entity;
62 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
63 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
66 static constexpr bool isSolDependent = (diffusionEnabled && diffusionIsSolDependent) ||
67 (heatConductionEnabled && heatConductionIsSolDependent);
71 : problemPtr_(&problem) {}
84 template<
class FluxVariablesCacheContainer,
class FluxVariablesCache,
class ElementVolumeVariables>
85 void fill(FluxVariablesCacheContainer& fluxVarsCacheContainer,
86 FluxVariablesCache& scvfFluxVarsCache,
87 const Element& element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const SubControlVolumeFace& scvf,
91 const bool forceUpdateAll =
false)
96 if constexpr (diffusionEnabled)
97 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
98 if constexpr (heatConductionEnabled)
99 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
103 if constexpr (diffusionEnabled && diffusionIsSolDependent)
104 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
105 if constexpr (heatConductionEnabled && heatConductionIsSolDependent)
106 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
112 const Problem& problem()
const
113 {
return *problemPtr_; }
117 template<
class FluxVariablesCache,
class ElementVolumeVariables>
118 void fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
119 const Element& element,
120 const FVElementGeometry& fvGeometry,
121 const ElementVolumeVariables& elemVolVars,
122 const SubControlVolumeFace& scvf)
124 using DiffusionType =
typename ElementVolumeVariables::VolumeVariables::MolecularDiffusionType;
125 using DiffusionFiller =
typename DiffusionType::Cache::Filler;
126 using FluidSystem =
typename ElementVolumeVariables::VolumeVariables::FluidSystem;
128 static constexpr int numPhases = ModelTraits::numFluidPhases();
129 static constexpr int numComponents = ModelTraits::numFluidComponents();
132 if constexpr (FluidSystem::isTracerFluidSystem())
133 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
134 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
135 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
137 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
138 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
139 if (compIdx != FluidSystem::getMainComponent(phaseIdx))
140 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
144 template<
class FluxVariablesCache,
class ElementVolumeVariables>
145 void fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
146 const Element& element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& elemVolVars,
149 const SubControlVolumeFace& scvf)
151 using HeatConductionType =
typename ElementVolumeVariables::VolumeVariables::HeatConductionType;
152 using HeatConductionFiller =
typename HeatConductionType::Cache::Filler;
155 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
158 const Problem* problemPtr_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
std::decay_t< decltype(std::declval< Problem >().gridGeometry())> GridGeometry
Definition: common/typetraits/problem.hh:45
Definition: scalarfluxvariablescachefiller.hh:39
FreeFlowScalarFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets the problem pointer.
Definition: scalarfluxvariablescachefiller.hh:70
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:85
Declares all properties used in Dumux.
Type traits for problem classes.