12#ifndef DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
13#define DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
26template<
class TypeTag,
class DiscretizationMethod>
35template<
class TypeTag>
39template<
class TypeTag>
45 using GridView =
typename GridGeometry::GridView;
46 using FVElementGeometry =
typename GridGeometry::LocalView;
47 using SubControlVolume =
typename GridGeometry::SubControlVolume;
48 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
51 using Element =
typename GridView::template Codim<0>::Entity;
53 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
54 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
55 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
57 static constexpr bool advectionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
58 static constexpr bool diffusionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
59 static constexpr bool heatConductionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
63 static constexpr bool isSolDependent = (advectionEnabled && advectionIsSolDependent) ||
64 (diffusionEnabled && diffusionIsSolDependent) ||
65 (heatConductionEnabled && heatConductionIsSolDependent);
69 : problemPtr_(&problem) {}
82 template<
class FluxVariablesCacheContainer,
class FluxVariablesCache>
83 void fill(FluxVariablesCacheContainer& fluxVarsCacheContainer,
84 FluxVariablesCache& scvfFluxVarsCache,
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const SubControlVolumeFace& scvf,
89 bool forceUpdateAll =
false)
94 if constexpr (advectionEnabled)
95 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
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 (advectionEnabled && advectionIsSolDependent)
104 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
105 if constexpr (diffusionEnabled && diffusionIsSolDependent)
106 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
107 if constexpr (heatConductionEnabled && heatConductionIsSolDependent)
108 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
114 const Problem& problem()
const
115 {
return *problemPtr_; }
118 template<
class FluxVariablesCache>
119 void fillAdvection_(FluxVariablesCache& scvfFluxVarsCache,
120 const Element& element,
121 const FVElementGeometry& fvGeometry,
122 const ElementVolumeVariables& elemVolVars,
123 const SubControlVolumeFace& scvf)
125 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
126 using AdvectionFiller =
typename AdvectionType::Cache::Filler;
129 AdvectionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
133 template<
class FluxVariablesCache>
134 void fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
135 const Element& element,
136 const FVElementGeometry& fvGeometry,
137 const ElementVolumeVariables& elemVolVars,
138 const SubControlVolumeFace& scvf)
140 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
141 using DiffusionFiller =
typename DiffusionType::Cache::Filler;
142 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
144 static constexpr int numPhases = ModelTraits::numFluidPhases();
145 static constexpr int numComponents = ModelTraits::numFluidComponents();
148 if constexpr (FluidSystem::isTracerFluidSystem())
149 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
150 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
151 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(),
element, fvGeometry, elemVolVars, scvf, *
this);
153 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
154 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
155 if (compIdx != FluidSystem::getMainComponent(phaseIdx))
156 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(),
element, fvGeometry, elemVolVars, scvf, *
this);
160 template<
class FluxVariablesCache>
161 void fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
162 const Element& element,
163 const FVElementGeometry& fvGeometry,
164 const ElementVolumeVariables& elemVolVars,
165 const SubControlVolumeFace& scvf)
167 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
168 using HeatConductionFiller =
typename HeatConductionType::Cache::Filler;
171 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
174 const Problem* problemPtr_;
178template<
class TypeTag>
184 using Element =
typename GridView::template Codim<0>::Entity;
188 using FVElementGeometry =
typename GridGeometry::LocalView;
189 using MpfaHelper =
typename GridGeometry::MpfaHelper;
190 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
196 using PrimaryDataHandle =
typename ElementFluxVariablesCache::PrimaryIvDataHandle;
197 using PrimaryLocalFaceData =
typename PrimaryInteractionVolume::Traits::LocalFaceData;
199 using SecondaryDataHandle =
typename ElementFluxVariablesCache::SecondaryIvDataHandle;
200 using SecondaryLocalFaceData =
typename SecondaryInteractionVolume::Traits::LocalFaceData;
202 static constexpr int dim = GridView::dimension;
203 static constexpr int dimWorld = GridView::dimensionworld;
205 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
206 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
207 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
209 static constexpr bool advectionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
210 static constexpr bool diffusionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
211 static constexpr bool heatConductionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
217 static constexpr bool isSolDependent =
true;
221 : problemPtr_(&problem) {}
234 template<
class FluxVarsCacheStorage,
class FluxVariablesCache,
class IVDataStorage>
235 void fill(FluxVarsCacheStorage& fluxVarsCacheStorage,
236 FluxVariablesCache& scvfFluxVarsCache,
237 IVDataStorage& ivDataStorage,
238 const FVElementGeometry& fvGeometry,
239 const ElementVolumeVariables& elemVolVars,
240 const SubControlVolumeFace& scvf,
241 bool forceUpdateAll =
false)
244 fvGeometryPtr_ = &fvGeometry;
245 elemVolVarsPtr_ = &elemVolVars;
246 const auto& gridGeometry = fvGeometry.gridGeometry();
251 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
256 const auto ivIndexInContainer = ivDataStorage.secondaryInteractionVolumes.size();
257 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().secondaryIndexSet(scvf);
258 ivDataStorage.secondaryInteractionVolumes.emplace_back();
259 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes.back();
260 secondaryIv_->bind(indexSet, problem(), fvGeometry);
263 ivDataStorage.secondaryDataHandles.emplace_back();
264 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles.back();
265 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
268 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
273 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
274 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes[ivIndexInContainer];
275 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles[ivIndexInContainer];
276 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
279 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
289 const auto ivIndexInContainer = ivDataStorage.primaryInteractionVolumes.size();
290 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().primaryIndexSet(scvf);
291 ivDataStorage.primaryInteractionVolumes.emplace_back();
292 primaryIv_ = &ivDataStorage.primaryInteractionVolumes.back();
293 primaryIv_->bind(indexSet, problem(), fvGeometry);
296 ivDataStorage.primaryDataHandles.emplace_back();
297 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles.back();
298 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
301 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
306 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
307 primaryIv_ = &ivDataStorage.primaryInteractionVolumes[ivIndexInContainer];
308 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles[ivIndexInContainer];
309 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
312 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
319 {
return *primaryIv_; }
323 {
return *secondaryIv_; }
327 {
return *primaryIvDataHandle_; }
331 {
return *secondaryIvDataHandle_; }
335 {
return *primaryLocalFaceData_; }
339 {
return *secondaryLocalFaceData_; }
343 const Problem& problem()
const {
return *problemPtr_; }
344 const FVElementGeometry& fvGeometry()
const {
return *fvGeometryPtr_; }
345 const ElementVolumeVariables& elemVolVars()
const {
return *elemVolVarsPtr_; }
348 template<
class FluxVariablesCache,
class FluxVarsCacheStorage,
class InteractionVolume>
349 void fillCachesInInteractionVolume_(FluxVarsCacheStorage& fluxVarsCacheStorage,
350 InteractionVolume& iv,
351 unsigned int ivIndexInContainer)
354 static constexpr bool isSecondary = MpfaHelper::considerSecondaryIVs()
355 && std::is_same_v<InteractionVolume, SecondaryInteractionVolume>;
360 const auto numGlobalScvfs = iv.localFaceData().size();
361 std::vector<const SubControlVolumeFace*> ivScvfs(numGlobalScvfs);
362 std::vector<FluxVariablesCache*> ivFluxVarCaches(numGlobalScvfs);
365 for (
const auto& d : iv.localFaceData())
368 const auto& scvfJ = fvGeometry().scvf(d.gridScvfIndex());
370 ivFluxVarCaches[i] = &fluxVarsCacheStorage[scvfJ];
371 ivFluxVarCaches[i]->setIvIndexInContainer(ivIndexInContainer);
372 ivFluxVarCaches[i]->setUpdateStatus(
true);
373 ivFluxVarCaches[i]->setSecondaryIvUsage(isSecondary);
374 ivFluxVarCaches[i]->setIvLocalFaceIndex(d.ivLocalScvfIndex());
376 if (d.isOutsideFace())
377 ivFluxVarCaches[i]->setIndexInOutsideFaces(d.scvfLocalOutsideScvfIndex());
381 if constexpr (advectionEnabled)
382 fillAdvection_(iv, ivScvfs, ivFluxVarCaches);
383 if constexpr (diffusionEnabled)
384 fillDiffusion_(iv, ivScvfs, ivFluxVarCaches);
385 if constexpr (heatConductionEnabled)
386 fillHeatConduction_(iv, ivScvfs, ivFluxVarCaches);
390 template<
class InteractionVolume,
class FluxVariablesCache>
391 void fillAdvection_(InteractionVolume& iv,
392 const std::vector<const SubControlVolumeFace*>& ivScvfs,
393 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
395 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
396 using AdvectionFiller =
typename AdvectionType::Cache::Filler;
399 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
403 if (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
406 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
407 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
409 else if (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
410 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
412 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
415 AdvectionFiller::fill(*ivFluxVarCaches[i],
417 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
426 template<
class InteractionVolume,
class FluxVariablesCache>
427 void fillDiffusion_(InteractionVolume& iv,
428 const std::vector<const SubControlVolumeFace*>& ivScvfs,
429 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
431 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
432 using DiffusionFiller =
typename DiffusionType::Cache::Filler;
434 static constexpr int numPhases = ModelTraits::numFluidPhases();
435 static constexpr int numComponents = ModelTraits::numFluidComponents();
437 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
439 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
441 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
442 if constexpr (!FluidSystem::isTracerFluidSystem())
443 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
447 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
451 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
454 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
455 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
457 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
458 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
460 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
463 DiffusionFiller::fill(*ivFluxVarCaches[i],
467 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
478 template<
class InteractionVolume,
class FluxVariablesCache>
479 void fillHeatConduction_(InteractionVolume& iv,
480 const std::vector<const SubControlVolumeFace*>& ivScvfs,
481 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
483 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
484 using HeatConductionFiller =
typename HeatConductionType::Cache::Filler;
487 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
491 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
494 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
495 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
497 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
498 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
500 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
503 HeatConductionFiller::fill(*ivFluxVarCaches[i],
505 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
514 template<
class InteractionVolume,
class DataHandle>
515 void prepareDataHandle_([[maybe_unused]] InteractionVolume& iv, [[maybe_unused]] DataHandle& handle, [[maybe_unused]]
bool forceUpdate)
518 if constexpr (advectionEnabled)
520 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
522 prepareAdvectionHandle_(iv, handle, forceUpdate);
526 if constexpr (diffusionEnabled)
528 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
530 prepareDiffusionHandles_(iv, handle, forceUpdate);
534 if constexpr (heatConductionEnabled)
536 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
538 prepareHeatConductionHandle_(iv, handle, forceUpdate);
543 template<
class InteractionVolume,
class DataHandle>
544 void prepareAdvectionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll)
547 using Traits =
typename InteractionVolume::Traits;
548 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
549 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
552 auto getK = [] (
const auto& volVars) {
return volVars.permeability(); };
555 if (forceUpdateAll || advectionIsSolDependent)
556 localAssembler.assembleMatrices(handle.advectionHandle(), iv, getK);
559 for (
unsigned int pIdx = 0; pIdx < ModelTraits::numFluidPhases(); ++pIdx)
562 handle.advectionHandle().setPhaseIndex(pIdx);
565 auto getRho = [pIdx] (
const auto& volVars) {
return volVars.density(pIdx); };
566 static const bool enableGravity = getParamFromGroup<bool>(problem().paramGroup(),
"Problem.EnableGravity");
568 localAssembler.assembleGravity(handle.advectionHandle(), iv, getRho);
571 auto getPressure = [pIdx] (
const auto& volVars) {
return volVars.pressure(pIdx); };
572 localAssembler.assembleU(handle.advectionHandle(), iv, getPressure);
577 template<
class InteractionVolume,
class DataHandle>
578 void prepareDiffusionHandles_(InteractionVolume& iv,
582 for (
unsigned int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
584 for (
unsigned int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
587 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
588 if constexpr (!FluidSystem::isTracerFluidSystem())
589 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
593 handle.diffusionHandle().setPhaseIndex(phaseIdx);
594 handle.diffusionHandle().setComponentIndex(compIdx);
596 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
599 using Traits =
typename InteractionVolume::Traits;
600 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
601 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
604 if (forceUpdateAll || diffusionIsSolDependent)
607 const auto getD = [&](
const auto& volVars)
609 if constexpr (FluidSystem::isTracerFluidSystem())
610 return volVars.effectiveDiffusionCoefficient(0, 0, compIdx);
612 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
617 static const auto zeroD = getParamFromGroup<Scalar>(
618 problem().paramGroup(),
619 "Mpfa.ZeroEffectiveDiffusionCoefficientThreshold",
625 const auto& scv = fvGeometry().scv(iv.localScv(0).gridScvIndex());
626 const auto& scvf = fvGeometry().scvf(iv.localScvf(0).gridScvfIndex());
627 const auto& vv = elemVolVars()[scv];
629 fvGeometry(), scvf, scv, zeroD, vv.extrusionFactor()
632 localAssembler.assembleMatrices(handle.diffusionHandle(), iv, getD, eps);
636 auto getMassOrMoleFraction = [phaseIdx, compIdx] (
const auto& volVars)
639 volVars.moleFraction(phaseIdx, compIdx);
642 localAssembler.assembleU(handle.diffusionHandle(), iv, getMassOrMoleFraction);
648 template<
class InteractionVolume,
class DataHandle>
649 void prepareHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll)
652 using Traits =
typename InteractionVolume::Traits;
653 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
654 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
657 auto getLambda = [] (
const auto& volVars) {
return volVars.effectiveThermalConductivity(); };
660 if (forceUpdateAll || heatConductionIsSolDependent)
661 localAssembler.assembleMatrices(handle.heatConductionHandle(), iv, getLambda);
664 auto getTemperature = [] (
const auto& volVars) {
return volVars.temperature(); };
665 localAssembler.assembleU(handle.heatConductionHandle(), iv, getTemperature);
668 const Problem* problemPtr_;
669 const FVElementGeometry* fvGeometryPtr_;
670 const ElementVolumeVariables* elemVolVarsPtr_;
675 PrimaryInteractionVolume* primaryIv_;
676 SecondaryInteractionVolume* secondaryIv_;
679 PrimaryDataHandle* primaryIvDataHandle_;
680 SecondaryDataHandle* secondaryIvDataHandle_;
686 const PrimaryLocalFaceData* primaryLocalFaceData_;
687 const SecondaryLocalFaceData* secondaryLocalFaceData_;
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets the problem pointer.
Definition: fluxvariablescachefiller.hh:68
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:83
const PrimaryInteractionVolume & primaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:318
const SecondaryDataHandle & secondaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:330
const PrimaryDataHandle & primaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:326
const SecondaryLocalFaceData & secondaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:338
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:235
const PrimaryLocalFaceData & primaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:334
const SecondaryInteractionVolume & secondaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:322
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets problem pointer.
Definition: fluxvariablescachefiller.hh:220
Definition: fluxvariablescachefiller.hh:27
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
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
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
constexpr CCMpfa ccmpfa
Definition: method.hh:146
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...