24#ifndef DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
25#define DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
38template<
class TypeTag,
class DiscretizationMethod>
47template<
class TypeTag>
51template<
class TypeTag>
57 using GridView =
typename GridGeometry::GridView;
58 using FVElementGeometry =
typename GridGeometry::LocalView;
59 using SubControlVolume =
typename GridGeometry::SubControlVolume;
60 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
63 using Element =
typename GridView::template Codim<0>::Entity;
65 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
66 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
67 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
69 static constexpr bool advectionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
70 static constexpr bool diffusionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
71 static constexpr bool heatConductionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
75 static constexpr bool isSolDependent = (advectionEnabled && advectionIsSolDependent) ||
76 (diffusionEnabled && diffusionIsSolDependent) ||
77 (heatConductionEnabled && heatConductionIsSolDependent);
81 : problemPtr_(&problem) {}
94 template<
class FluxVariablesCacheContainer,
class FluxVariablesCache>
95 void fill(FluxVariablesCacheContainer& fluxVarsCacheContainer,
96 FluxVariablesCache& scvfFluxVarsCache,
97 const Element& element,
98 const FVElementGeometry& fvGeometry,
99 const ElementVolumeVariables& elemVolVars,
100 const SubControlVolumeFace& scvf,
101 bool forceUpdateAll =
false)
106 if constexpr (advectionEnabled)
107 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
108 if constexpr (diffusionEnabled)
109 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
110 if constexpr (heatConductionEnabled)
111 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
115 if constexpr (advectionEnabled && advectionIsSolDependent)
116 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
117 if constexpr (diffusionEnabled && diffusionIsSolDependent)
118 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
119 if constexpr (heatConductionEnabled && heatConductionIsSolDependent)
120 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
126 const Problem& problem()
const
127 {
return *problemPtr_; }
130 template<
class FluxVariablesCache>
131 void fillAdvection_(FluxVariablesCache& scvfFluxVarsCache,
132 const Element& element,
133 const FVElementGeometry& fvGeometry,
134 const ElementVolumeVariables& elemVolVars,
135 const SubControlVolumeFace& scvf)
137 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
138 using AdvectionFiller =
typename AdvectionType::Cache::Filler;
141 AdvectionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
145 template<
class FluxVariablesCache>
146 void fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
147 const Element& element,
148 const FVElementGeometry& fvGeometry,
149 const ElementVolumeVariables& elemVolVars,
150 const SubControlVolumeFace& scvf)
152 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
153 using DiffusionFiller =
typename DiffusionType::Cache::Filler;
154 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
156 static constexpr int numPhases = ModelTraits::numFluidPhases();
157 static constexpr int numComponents = ModelTraits::numFluidComponents();
160 if constexpr (FluidSystem::isTracerFluidSystem())
161 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
162 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
163 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(),
element, fvGeometry, elemVolVars, scvf, *
this);
165 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
166 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
167 if (compIdx != FluidSystem::getMainComponent(phaseIdx))
168 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(),
element, fvGeometry, elemVolVars, scvf, *
this);
172 template<
class FluxVariablesCache>
173 void fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
174 const Element& element,
175 const FVElementGeometry& fvGeometry,
176 const ElementVolumeVariables& elemVolVars,
177 const SubControlVolumeFace& scvf)
179 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
180 using HeatConductionFiller =
typename HeatConductionType::Cache::Filler;
183 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
186 const Problem* problemPtr_;
190template<
class TypeTag>
196 using Element =
typename GridView::template Codim<0>::Entity;
200 using FVElementGeometry =
typename GridGeometry::LocalView;
201 using MpfaHelper =
typename GridGeometry::MpfaHelper;
202 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
208 using PrimaryDataHandle =
typename ElementFluxVariablesCache::PrimaryIvDataHandle;
209 using PrimaryLocalFaceData =
typename PrimaryInteractionVolume::Traits::LocalFaceData;
211 using SecondaryDataHandle =
typename ElementFluxVariablesCache::SecondaryIvDataHandle;
212 using SecondaryLocalFaceData =
typename SecondaryInteractionVolume::Traits::LocalFaceData;
214 static constexpr int dim = GridView::dimension;
215 static constexpr int dimWorld = GridView::dimensionworld;
217 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
218 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
219 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
221 static constexpr bool advectionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
222 static constexpr bool diffusionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
223 static constexpr bool heatConductionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
229 static constexpr bool isSolDependent =
true;
233 : problemPtr_(&problem) {}
246 template<
class FluxVarsCacheStorage,
class FluxVariablesCache,
class IVDataStorage>
247 void fill(FluxVarsCacheStorage& fluxVarsCacheStorage,
248 FluxVariablesCache& scvfFluxVarsCache,
249 IVDataStorage& ivDataStorage,
250 const FVElementGeometry& fvGeometry,
251 const ElementVolumeVariables& elemVolVars,
252 const SubControlVolumeFace& scvf,
253 bool forceUpdateAll =
false)
256 fvGeometryPtr_ = &fvGeometry;
257 elemVolVarsPtr_ = &elemVolVars;
258 const auto& gridGeometry = fvGeometry.gridGeometry();
263 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
268 const auto ivIndexInContainer = ivDataStorage.secondaryInteractionVolumes.size();
269 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().secondaryIndexSet(scvf);
270 ivDataStorage.secondaryInteractionVolumes.emplace_back();
271 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes.back();
272 secondaryIv_->bind(indexSet, problem(), fvGeometry);
275 ivDataStorage.secondaryDataHandles.emplace_back();
276 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles.back();
277 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
280 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
285 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
286 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes[ivIndexInContainer];
287 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles[ivIndexInContainer];
288 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
291 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
301 const auto ivIndexInContainer = ivDataStorage.primaryInteractionVolumes.size();
302 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().primaryIndexSet(scvf);
303 ivDataStorage.primaryInteractionVolumes.emplace_back();
304 primaryIv_ = &ivDataStorage.primaryInteractionVolumes.back();
305 primaryIv_->bind(indexSet, problem(), fvGeometry);
308 ivDataStorage.primaryDataHandles.emplace_back();
309 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles.back();
310 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
313 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
318 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
319 primaryIv_ = &ivDataStorage.primaryInteractionVolumes[ivIndexInContainer];
320 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles[ivIndexInContainer];
321 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
324 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
331 {
return *primaryIv_; }
335 {
return *secondaryIv_; }
339 {
return *primaryIvDataHandle_; }
343 {
return *secondaryIvDataHandle_; }
347 {
return *primaryLocalFaceData_; }
351 {
return *secondaryLocalFaceData_; }
355 const Problem& problem()
const {
return *problemPtr_; }
356 const FVElementGeometry& fvGeometry()
const {
return *fvGeometryPtr_; }
357 const ElementVolumeVariables& elemVolVars()
const {
return *elemVolVarsPtr_; }
360 template<
class FluxVariablesCache,
class FluxVarsCacheStorage,
class InteractionVolume>
361 void fillCachesInInteractionVolume_(FluxVarsCacheStorage& fluxVarsCacheStorage,
362 InteractionVolume& iv,
363 unsigned int ivIndexInContainer)
366 static constexpr bool isSecondary = MpfaHelper::considerSecondaryIVs()
367 && std::is_same_v<InteractionVolume, SecondaryInteractionVolume>;
372 const auto numGlobalScvfs = iv.localFaceData().size();
373 std::vector<const SubControlVolumeFace*> ivScvfs(numGlobalScvfs);
374 std::vector<FluxVariablesCache*> ivFluxVarCaches(numGlobalScvfs);
377 for (
const auto& d : iv.localFaceData())
380 const auto& scvfJ = fvGeometry().scvf(d.gridScvfIndex());
382 ivFluxVarCaches[i] = &fluxVarsCacheStorage[scvfJ];
383 ivFluxVarCaches[i]->setIvIndexInContainer(ivIndexInContainer);
384 ivFluxVarCaches[i]->setUpdateStatus(
true);
385 ivFluxVarCaches[i]->setSecondaryIvUsage(isSecondary);
386 ivFluxVarCaches[i]->setIvLocalFaceIndex(d.ivLocalScvfIndex());
388 if (d.isOutsideFace())
389 ivFluxVarCaches[i]->setIndexInOutsideFaces(d.scvfLocalOutsideScvfIndex());
393 if constexpr (advectionEnabled)
394 fillAdvection_(iv, ivScvfs, ivFluxVarCaches);
395 if constexpr (diffusionEnabled)
396 fillDiffusion_(iv, ivScvfs, ivFluxVarCaches);
397 if constexpr (heatConductionEnabled)
398 fillHeatConduction_(iv, ivScvfs, ivFluxVarCaches);
402 template<
class InteractionVolume,
class FluxVariablesCache>
403 void fillAdvection_(InteractionVolume& iv,
404 const std::vector<const SubControlVolumeFace*>& ivScvfs,
405 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
407 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
408 using AdvectionFiller =
typename AdvectionType::Cache::Filler;
411 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
415 if (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
418 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
419 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
421 else if (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
422 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
424 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
427 AdvectionFiller::fill(*ivFluxVarCaches[i],
429 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
438 template<
class InteractionVolume,
class FluxVariablesCache>
439 void fillDiffusion_(InteractionVolume& iv,
440 const std::vector<const SubControlVolumeFace*>& ivScvfs,
441 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
443 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
444 using DiffusionFiller =
typename DiffusionType::Cache::Filler;
446 static constexpr int numPhases = ModelTraits::numFluidPhases();
447 static constexpr int numComponents = ModelTraits::numFluidComponents();
449 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
451 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
453 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
454 if constexpr (!FluidSystem::isTracerFluidSystem())
455 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
459 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
463 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
466 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
467 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
469 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
470 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
472 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
475 DiffusionFiller::fill(*ivFluxVarCaches[i],
479 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
490 template<
class InteractionVolume,
class FluxVariablesCache>
491 void fillHeatConduction_(InteractionVolume& iv,
492 const std::vector<const SubControlVolumeFace*>& ivScvfs,
493 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
495 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
496 using HeatConductionFiller =
typename HeatConductionType::Cache::Filler;
499 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
503 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
506 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
507 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
509 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
510 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
512 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
515 HeatConductionFiller::fill(*ivFluxVarCaches[i],
517 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
526 template<
class InteractionVolume,
class DataHandle>
527 void prepareDataHandle_([[maybe_unused]] InteractionVolume& iv, [[maybe_unused]] DataHandle& handle, [[maybe_unused]]
bool forceUpdate)
530 if constexpr (advectionEnabled)
532 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
534 prepareAdvectionHandle_(iv, handle, forceUpdate);
538 if constexpr (diffusionEnabled)
540 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
542 prepareDiffusionHandles_(iv, handle, forceUpdate);
546 if constexpr (heatConductionEnabled)
548 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
550 prepareHeatConductionHandle_(iv, handle, forceUpdate);
555 template<
class InteractionVolume,
class DataHandle>
556 void prepareAdvectionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll)
559 using Traits =
typename InteractionVolume::Traits;
560 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
561 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
564 auto getK = [] (
const auto& volVars) {
return volVars.permeability(); };
567 if (forceUpdateAll || advectionIsSolDependent)
568 localAssembler.assembleMatrices(handle.advectionHandle(), iv, getK);
571 for (
unsigned int pIdx = 0; pIdx < ModelTraits::numFluidPhases(); ++pIdx)
574 handle.advectionHandle().setPhaseIndex(pIdx);
577 auto getRho = [pIdx] (
const auto& volVars) {
return volVars.density(pIdx); };
578 static const bool enableGravity = getParamFromGroup<bool>(problem().paramGroup(),
"Problem.EnableGravity");
580 localAssembler.assembleGravity(handle.advectionHandle(), iv, getRho);
583 auto getPressure = [pIdx] (
const auto& volVars) {
return volVars.pressure(pIdx); };
584 localAssembler.assembleU(handle.advectionHandle(), iv, getPressure);
589 template<
class InteractionVolume,
class DataHandle>
590 void prepareDiffusionHandles_(InteractionVolume& iv,
594 for (
unsigned int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
596 for (
unsigned int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
599 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
600 if constexpr (!FluidSystem::isTracerFluidSystem())
601 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
605 handle.diffusionHandle().setPhaseIndex(phaseIdx);
606 handle.diffusionHandle().setComponentIndex(compIdx);
608 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
611 using Traits =
typename InteractionVolume::Traits;
612 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
613 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
616 if (forceUpdateAll || diffusionIsSolDependent)
619 const auto getD = [&](
const auto& volVars)
621 if constexpr (FluidSystem::isTracerFluidSystem())
622 return volVars.effectiveDiffusionCoefficient(0, 0, compIdx);
624 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
629 static const auto zeroD = getParamFromGroup<Scalar>(
630 problem().paramGroup(),
631 "Mpfa.ZeroEffectiveDiffusionCoefficientThreshold",
637 const auto& scv = fvGeometry().scv(iv.localScv(0).gridScvIndex());
638 const auto& scvf = fvGeometry().scvf(iv.localScvf(0).gridScvfIndex());
639 const auto& vv = elemVolVars()[scv];
641 scvf, scv, zeroD, vv.extrusionFactor()
644 localAssembler.assembleMatrices(handle.diffusionHandle(), iv, getD, eps);
648 auto getMassOrMoleFraction = [phaseIdx, compIdx] (
const auto& volVars)
651 volVars.moleFraction(phaseIdx, compIdx);
654 localAssembler.assembleU(handle.diffusionHandle(), iv, getMassOrMoleFraction);
660 template<
class InteractionVolume,
class DataHandle>
661 void prepareHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll)
664 using Traits =
typename InteractionVolume::Traits;
665 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
666 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
669 auto getLambda = [] (
const auto& volVars) {
return volVars.effectiveThermalConductivity(); };
672 if (forceUpdateAll || heatConductionIsSolDependent)
673 localAssembler.assembleMatrices(handle.heatConductionHandle(), iv, getLambda);
676 auto getTemperature = [] (
const auto& volVars) {
return volVars.temperature(); };
677 localAssembler.assembleU(handle.heatConductionHandle(), iv, getTemperature);
680 const Problem* problemPtr_;
681 const FVElementGeometry* fvGeometryPtr_;
682 const ElementVolumeVariables* elemVolVarsPtr_;
687 PrimaryInteractionVolume* primaryIv_;
688 SecondaryInteractionVolume* secondaryIv_;
691 PrimaryDataHandle* primaryIvDataHandle_;
692 SecondaryDataHandle* secondaryIvDataHandle_;
698 const PrimaryLocalFaceData* primaryLocalFaceData_;
699 const SecondaryLocalFaceData* secondaryLocalFaceData_;
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...
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename 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:47
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
constexpr CCMpfa ccmpfa
Definition: method.hh:138
Definition: fluxvariablescachefiller.hh:39
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets the problem pointer.
Definition: fluxvariablescachefiller.hh:80
void fill(FluxVariablesCacheContainer &fluxVarsCacheContainer, FluxVariablesCache &scvfFluxVarsCache, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, bool forceUpdateAll=false)
function to fill the flux variables caches
Definition: fluxvariablescachefiller.hh:95
const PrimaryInteractionVolume & primaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:330
const SecondaryDataHandle & secondaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:342
const PrimaryDataHandle & primaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:338
const SecondaryLocalFaceData & secondaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:350
void fill(FluxVarsCacheStorage &fluxVarsCacheStorage, FluxVariablesCache &scvfFluxVarsCache, IVDataStorage &ivDataStorage, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, bool forceUpdateAll=false)
function to fill the flux variables caches
Definition: fluxvariablescachefiller.hh:247
const PrimaryLocalFaceData & primaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:346
const SecondaryInteractionVolume & secondaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:334
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets problem pointer.
Definition: fluxvariablescachefiller.hh:232
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...