24#ifndef DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
25#define DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
37template<
class TypeTag, DiscretizationMethod discMethod>
46template<
class TypeTag>
50template<
class TypeTag>
57 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
61 using Element =
typename GridView::template Codim<0>::Entity;
63 static constexpr bool doAdvection = ModelTraits::enableAdvection();
64 static constexpr bool doDiffusion = ModelTraits::enableMolecularDiffusion();
65 static constexpr bool doHeatConduction = ModelTraits::enableEnergyBalance();
67 static constexpr bool soldependentAdvection = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
68 static constexpr bool soldependentDiffusion = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
69 static constexpr bool soldependentHeatConduction = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
73 static constexpr bool isSolDependent = (doAdvection && soldependentAdvection) ||
74 (doDiffusion && soldependentDiffusion) ||
75 (doHeatConduction && soldependentHeatConduction);
79 : problemPtr_(&problem) {}
92 template<
class FluxVariablesCacheContainer,
class FluxVariablesCache>
93 void fill(FluxVariablesCacheContainer& fluxVarsCacheContainer,
94 FluxVariablesCache& scvfFluxVarsCache,
95 const Element& element,
96 const FVElementGeometry& fvGeometry,
97 const ElementVolumeVariables& elemVolVars,
98 const SubControlVolumeFace& scvf,
99 bool forceUpdateAll =
false)
104 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
105 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
106 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
110 if (doAdvection && soldependentAdvection)
111 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
112 if (doDiffusion && soldependentDiffusion)
113 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
114 if (doHeatConduction && soldependentHeatConduction)
115 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
121 const Problem& problem()
const
122 {
return *problemPtr_; }
125 template<
class FluxVariablesCache,
bool advectionEnabled = doAdvection>
126 typename std::enable_if<advectionEnabled>::type
127 fillAdvection_(FluxVariablesCache& scvfFluxVarsCache,
128 const Element& element,
129 const FVElementGeometry& fvGeometry,
130 const ElementVolumeVariables& elemVolVars,
131 const SubControlVolumeFace& scvf)
133 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
134 using AdvectionFiller =
typename AdvectionType::Cache::Filler;
137 AdvectionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
141 template<
class FluxVariablesCache,
bool advectionEnabled = doAdvection>
142 typename std::enable_if<!advectionEnabled>::type
143 fillAdvection_(FluxVariablesCache& scvfFluxVarsCache,
144 const Element& element,
145 const FVElementGeometry& fvGeometry,
146 const ElementVolumeVariables& elemVolVars,
147 const SubControlVolumeFace& scvf)
151 template<
class FluxVariablesCache,
bool diffusionEnabled = doDiffusion>
152 typename std::enable_if<diffusionEnabled>::type
153 fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
154 const Element& element,
155 const FVElementGeometry& fvGeometry,
156 const ElementVolumeVariables& elemVolVars,
157 const SubControlVolumeFace& scvf)
159 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
160 using DiffusionFiller =
typename DiffusionType::Cache::Filler;
161 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
163 static constexpr int numPhases = ModelTraits::numFluidPhases();
164 static constexpr int numComponents = ModelTraits::numFluidComponents();
167 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
168 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
169 if (compIdx != FluidSystem::getMainComponent(phaseIdx))
170 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
174 template<
class FluxVariablesCache,
bool diffusionEnabled = doDiffusion>
175 typename std::enable_if<!diffusionEnabled>::type
176 fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
177 const Element& element,
178 const FVElementGeometry& fvGeometry,
179 const ElementVolumeVariables& elemVolVars,
180 const SubControlVolumeFace& scvf)
184 template<
class FluxVariablesCache,
bool heatConductionEnabled = doHeatConduction>
185 typename std::enable_if<heatConductionEnabled>::type
186 fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
187 const Element& element,
188 const FVElementGeometry& fvGeometry,
189 const ElementVolumeVariables& elemVolVars,
190 const SubControlVolumeFace& scvf)
192 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
193 using HeatConductionFiller =
typename HeatConductionType::Cache::Filler;
196 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *
this);
200 template<
class FluxVariablesCache,
bool heatConductionEnabled = doHeatConduction>
201 typename std::enable_if<!heatConductionEnabled>::type
202 fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
203 const Element& element,
204 const FVElementGeometry& fvGeometry,
205 const ElementVolumeVariables& elemVolVars,
206 const SubControlVolumeFace& scvf)
209 const Problem* problemPtr_;
213template<
class TypeTag>
219 using Element =
typename GridView::template Codim<0>::Entity;
222 using FVElementGeometry =
typename GridGeometry::LocalView;
223 using MpfaHelper =
typename GridGeometry::MpfaHelper;
224 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
229 using PrimaryDataHandle =
typename ElementFluxVariablesCache::PrimaryIvDataHandle;
230 using PrimaryLocalFaceData =
typename PrimaryInteractionVolume::Traits::LocalFaceData;
232 using SecondaryDataHandle =
typename ElementFluxVariablesCache::SecondaryIvDataHandle;
233 using SecondaryLocalFaceData =
typename SecondaryInteractionVolume::Traits::LocalFaceData;
235 static constexpr int dim = GridView::dimension;
236 static constexpr int dimWorld = GridView::dimensionworld;
238 static constexpr bool doAdvection = ModelTraits::enableAdvection();
239 static constexpr bool doDiffusion = ModelTraits::enableMolecularDiffusion();
240 static constexpr bool doHeatConduction = ModelTraits::enableEnergyBalance();
242 static constexpr bool soldependentAdvection = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
243 static constexpr bool soldependentDiffusion = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
244 static constexpr bool soldependentHeatConduction = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
250 static constexpr bool isSolDependent =
true;
254 : problemPtr_(&problem) {}
268 template<
class FluxVarsCacheStorage,
class FluxVariablesCache,
class IVDataStorage>
269 void fill(FluxVarsCacheStorage& fluxVarsCacheStorage,
270 FluxVariablesCache& scvfFluxVarsCache,
271 IVDataStorage& ivDataStorage,
272 const Element& element,
273 const FVElementGeometry& fvGeometry,
274 const ElementVolumeVariables& elemVolVars,
275 const SubControlVolumeFace& scvf,
276 bool forceUpdateAll =
false)
279 elementPtr_ = &element;
280 fvGeometryPtr_ = &fvGeometry;
281 elemVolVarsPtr_ = &elemVolVars;
284 const auto& gridGeometry = fvGeometry.gridGeometry();
285 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
290 const auto ivIndexInContainer = ivDataStorage.secondaryInteractionVolumes.size();
293 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().secondaryIndexSet(scvf);
294 ivDataStorage.secondaryInteractionVolumes.emplace_back();
295 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes.back();
296 secondaryIv_->bind(indexSet, problem(), fvGeometry);
299 ivDataStorage.secondaryDataHandles.emplace_back();
300 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles.back();
303 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, *secondaryIvDataHandle_, ivIndexInContainer,
true);
307 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
308 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes[ivIndexInContainer];
309 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles[ivIndexInContainer];
312 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, *secondaryIvDataHandle_, ivIndexInContainer);
320 const auto ivIndexInContainer = ivDataStorage.primaryInteractionVolumes.size();
323 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().primaryIndexSet(scvf);
324 ivDataStorage.primaryInteractionVolumes.emplace_back();
325 primaryIv_ = &ivDataStorage.primaryInteractionVolumes.back();
326 primaryIv_->bind(indexSet, problem(), fvGeometry);
329 ivDataStorage.primaryDataHandles.emplace_back();
330 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles.back();
333 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, *primaryIvDataHandle_, ivIndexInContainer,
true);
337 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
338 primaryIv_ = &ivDataStorage.primaryInteractionVolumes[ivIndexInContainer];
339 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles[ivIndexInContainer];
342 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, *primaryIvDataHandle_, ivIndexInContainer);
349 {
return *primaryIv_; }
353 {
return *secondaryIv_; }
357 {
return *primaryIvDataHandle_; }
361 {
return *secondaryIvDataHandle_; }
365 {
return *primaryLocalFaceData_; }
369 {
return *secondaryLocalFaceData_; }
373 const Problem& problem()
const {
return *problemPtr_; }
374 const Element& element()
const {
return *elementPtr_; }
375 const FVElementGeometry& fvGeometry()
const {
return *fvGeometryPtr_; }
376 const ElementVolumeVariables& elemVolVars()
const {
return *elemVolVarsPtr_; }
379 template<
class FluxVariablesCache,
class FluxVarsCacheStorage,
class InteractionVolume,
class DataHandle>
380 void fillCachesInInteractionVolume_(FluxVarsCacheStorage& fluxVarsCacheStorage,
381 InteractionVolume& iv,
383 unsigned int ivIndexInContainer,
384 bool forceUpdateAll =
false)
387 static constexpr bool isSecondary = MpfaHelper::considerSecondaryIVs()
388 && std::is_same<InteractionVolume, SecondaryInteractionVolume>::value;
393 const auto numGlobalScvfs = iv.localFaceData().size();
394 std::vector<const SubControlVolumeFace*> ivScvfs(numGlobalScvfs);
395 std::vector<FluxVariablesCache*> ivFluxVarCaches(numGlobalScvfs);
398 for (
const auto& d : iv.localFaceData())
401 const auto& scvfJ = fvGeometry().scvf(d.gridScvfIndex());
403 ivFluxVarCaches[i] = &fluxVarsCacheStorage[scvfJ];
404 ivFluxVarCaches[i]->setIvIndexInContainer(ivIndexInContainer);
405 ivFluxVarCaches[i]->setUpdateStatus(
true);
406 ivFluxVarCaches[i]->setSecondaryIvUsage(isSecondary);
407 ivFluxVarCaches[i]->setIvLocalFaceIndex(d.ivLocalScvfIndex());
409 if (d.isOutsideFace())
410 ivFluxVarCaches[i]->setIndexInOutsideFaces(d.scvfLocalOutsideScvfIndex());
414 fillAdvection_(iv, handle, ivScvfs, ivFluxVarCaches, forceUpdateAll);
415 fillDiffusion_(iv, handle, ivScvfs, ivFluxVarCaches, forceUpdateAll);
416 fillHeatConduction_(iv, handle, ivScvfs, ivFluxVarCaches, forceUpdateAll);
420 template<
class InteractionVolume,
422 class FluxVariablesCache,
423 bool enableAdvection = doAdvection,
424 typename std::enable_if_t<enableAdvection, int> = 0 >
425 void fillAdvection_(InteractionVolume& iv,
427 const std::vector<const SubControlVolumeFace*>& ivScvfs,
428 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
429 bool forceUpdateAll =
false)
431 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
432 using AdvectionFiller =
typename AdvectionType::Cache::Filler;
435 fillAdvectionHandle_(iv, handle, forceUpdateAll);
438 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
442 if (std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value)
445 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
446 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
448 else if (std::is_same<InteractionVolume, PrimaryInteractionVolume>::value)
449 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
451 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
454 AdvectionFiller::fill(*ivFluxVarCaches[i],
456 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
465 template<
class InteractionVolume,
467 class FluxVariablesCache,
468 bool enableAdvection = doAdvection,
469 typename std::enable_if_t<!enableAdvection, int> = 0 >
470 void fillAdvection_(InteractionVolume& iv,
472 const std::vector<const SubControlVolumeFace*>& ivScvfs,
473 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
474 bool forceUpdateAll =
false)
478 template<
class InteractionVolume,
480 class FluxVariablesCache,
481 bool enableDiffusion = doDiffusion,
482 typename std::enable_if_t<enableDiffusion, int> = 0 >
483 void fillDiffusion_(InteractionVolume& iv,
485 const std::vector<const SubControlVolumeFace*>& ivScvfs,
486 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
487 bool forceUpdateAll =
false)
489 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
490 using DiffusionFiller =
typename DiffusionType::Cache::Filler;
492 static constexpr int numPhases = ModelTraits::numFluidPhases();
493 static constexpr int numComponents = ModelTraits::numFluidComponents();
495 for (
unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
497 handle.diffusionHandle().setPhaseIndex(phaseIdx);
498 for (
unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
500 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
501 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
505 handle.diffusionHandle().setComponentIndex(compIdx);
506 fillDiffusionHandle_(iv, handle, forceUpdateAll, phaseIdx, compIdx);
509 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
513 if (std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value)
516 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
517 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
519 else if (std::is_same<InteractionVolume, PrimaryInteractionVolume>::value)
520 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
522 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
525 DiffusionFiller::fill(*ivFluxVarCaches[i],
529 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
540 template<
class InteractionVolume,
542 class FluxVariablesCache,
543 bool enableDiffusion = doDiffusion,
544 typename std::enable_if_t<!enableDiffusion, int> = 0 >
545 void fillDiffusion_(InteractionVolume& iv,
547 const std::vector<const SubControlVolumeFace*>& ivScvfs,
548 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
549 bool forceUpdateAll =
false)
553 template<
class InteractionVolume,
555 class FluxVariablesCache,
556 bool enableHeatConduction = doHeatConduction,
557 typename std::enable_if_t<enableHeatConduction, int> = 0 >
558 void fillHeatConduction_(InteractionVolume& iv,
560 const std::vector<const SubControlVolumeFace*>& ivScvfs,
561 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
562 bool forceUpdateAll =
false)
564 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
565 using HeatConductionFiller =
typename HeatConductionType::Cache::Filler;
568 fillHeatConductionHandle_(iv, handle, forceUpdateAll);
571 for (
unsigned int i = 0; i < iv.localFaceData().size(); ++i)
575 if (std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value)
578 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
579 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
581 else if (std::is_same<InteractionVolume, PrimaryInteractionVolume>::value)
582 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
584 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
587 HeatConductionFiller::fill(*ivFluxVarCaches[i],
589 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
598 template<
class InteractionVolume,
600 class FluxVariablesCache,
601 bool enableHeatConduction = doHeatConduction,
602 typename std::enable_if_t<!enableHeatConduction, int> = 0 >
603 void fillHeatConduction_(InteractionVolume& iv,
605 const std::vector<const SubControlVolumeFace*>& ivScvfs,
606 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
607 bool forceUpdateAll =
false)
611 template<
class InteractionVolume,
613 class AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>,
614 typename std::enable_if_t<AdvectionType::discMethod == DiscretizationMethod::ccmpfa, int> = 0 >
615 void fillAdvectionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll)
617 using LambdaFactory = TensorLambdaFactory<DiscretizationMethod::ccmpfa>;
620 using Traits =
typename InteractionVolume::Traits;
621 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
622 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
625 if (forceUpdateAll || soldependentAdvection)
626 localAssembler.assembleMatrices(handle.advectionHandle(), iv, LambdaFactory::getAdvectionLambda());
629 for (
unsigned int pIdx = 0; pIdx < ModelTraits::numFluidPhases(); ++pIdx)
632 handle.advectionHandle().setPhaseIndex(pIdx);
635 auto getRho = [pIdx] (
const auto& volVars) {
return volVars.density(pIdx); };
636 static const bool enableGravity = getParamFromGroup<bool>(problem().paramGroup(),
"Problem.EnableGravity");
638 localAssembler.assembleGravity(handle.advectionHandle(), iv, getRho);
641 auto getPressure = [pIdx] (
const auto& volVars) {
return volVars.pressure(pIdx); };
642 localAssembler.assembleU(handle.advectionHandle(), iv, getPressure);
647 template<
class InteractionVolume,
649 class DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>,
650 typename std::enable_if_t<DiffusionType::discMethod == DiscretizationMethod::ccmpfa, int> = 0 >
651 void fillDiffusionHandle_(InteractionVolume& iv,
654 int phaseIdx,
int compIdx)
656 using LambdaFactory = TensorLambdaFactory<DiscretizationMethod::ccmpfa>;
659 using Traits =
typename InteractionVolume::Traits;
660 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
661 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
664 if (forceUpdateAll || soldependentDiffusion)
665 localAssembler.assembleMatrices(handle.diffusionHandle(),
667 LambdaFactory::getDiffusionLambda(phaseIdx, compIdx));
670 auto getMassOrMoleFraction = [phaseIdx, compIdx] (
const auto& volVars) {
return (DiffusionType::referenceSystemFormulation() ==
ReferenceSystemFormulation::massAveraged) ?volVars.massFraction(phaseIdx, compIdx) : volVars.moleFraction(phaseIdx, compIdx); };
671 localAssembler.assembleU(handle.diffusionHandle(), iv, getMassOrMoleFraction);
675 template<
class InteractionVolume,
677 class HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>,
678 typename std::enable_if_t<HeatConductionType::discMethod == DiscretizationMethod::ccmpfa, int> = 0 >
679 void fillHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll)
681 using LambdaFactory = TensorLambdaFactory<DiscretizationMethod::ccmpfa>;
682 using ThermCondModel = GetPropType<TypeTag, Properties::ThermalConductivityModel>;
685 using Traits =
typename InteractionVolume::Traits;
686 using IvLocalAssembler =
typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
687 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
690 if (forceUpdateAll || soldependentHeatConduction)
691 localAssembler.assembleMatrices(handle.heatConductionHandle(),
693 LambdaFactory::template getHeatConductionLambda<ThermCondModel>());
696 auto getTemperature = [] (
const auto& volVars) {
return volVars.temperature(); };
697 localAssembler.assembleU(handle.heatConductionHandle(), iv, getTemperature);
701 template<
class InteractionVolume,
703 class AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>,
704 typename std::enable_if_t<AdvectionType::discMethod != DiscretizationMethod::ccmpfa, int> = 0 >
705 void fillAdvectionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll) {}
708 template<
class InteractionVolume,
710 class DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>,
711 typename std::enable_if_t<DiffusionType::discMethod != DiscretizationMethod::ccmpfa, int> = 0 >
712 void fillDiffusionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll,
int phaseIdx,
int compIdx) {}
715 template<
class InteractionVolume,
717 class HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>,
718 typename std::enable_if_t<HeatConductionType::discMethod != DiscretizationMethod::ccmpfa, int> = 0 >
719 void fillHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle,
bool forceUpdateAll) {}
721 const Problem* problemPtr_;
722 const Element* elementPtr_;
723 const FVElementGeometry* fvGeometryPtr_;
724 const ElementVolumeVariables* elemVolVarsPtr_;
729 PrimaryInteractionVolume* primaryIv_;
730 SecondaryInteractionVolume* secondaryIv_;
733 PrimaryDataHandle* primaryIvDataHandle_;
734 SecondaryDataHandle* secondaryIvDataHandle_;
740 const PrimaryLocalFaceData* primaryLocalFaceData_;
741 const SecondaryLocalFaceData* secondaryLocalFaceData_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper class to be used to obtain lambda functions for the tensors involved in the laws that describe...
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Definition: porousmediumflow/fluxvariablescachefiller.hh:38
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets the problem pointer.
Definition: porousmediumflow/fluxvariablescachefiller.hh:78
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: porousmediumflow/fluxvariablescachefiller.hh:93
const PrimaryLocalFaceData & primaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: porousmediumflow/fluxvariablescachefiller.hh:364
void fill(FluxVarsCacheStorage &fluxVarsCacheStorage, FluxVariablesCache &scvfFluxVarsCache, IVDataStorage &ivDataStorage, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, bool forceUpdateAll=false)
function to fill the flux variables caches
Definition: porousmediumflow/fluxvariablescachefiller.hh:269
const PrimaryInteractionVolume & primaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: porousmediumflow/fluxvariablescachefiller.hh:348
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets problem pointer.
Definition: porousmediumflow/fluxvariablescachefiller.hh:253
const PrimaryDataHandle & primaryIvDataHandle() const
returns the stored data handle pointer
Definition: porousmediumflow/fluxvariablescachefiller.hh:356
const SecondaryDataHandle & secondaryIvDataHandle() const
returns the stored data handle pointer
Definition: porousmediumflow/fluxvariablescachefiller.hh:360
const SecondaryLocalFaceData & secondaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: porousmediumflow/fluxvariablescachefiller.hh:368
const SecondaryInteractionVolume & secondaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: porousmediumflow/fluxvariablescachefiller.hh:352
Declares all properties used in Dumux.