version 3.8
fluxvariablescachefiller.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
13#define DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
14
17
22
23namespace Dumux {
24
25// forward declaration
26template<class TypeTag, class DiscretizationMethod>
28
35template<class TypeTag>
37
39template<class TypeTag>
40class PorousMediumFluxVariablesCacheFillerImplementation<TypeTag, DiscretizationMethods::CCTpfa>
41{
45 using GridView = typename GridGeometry::GridView;
46 using FVElementGeometry = typename GridGeometry::LocalView;
47 using SubControlVolume = typename GridGeometry::SubControlVolume;
48 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
49 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
50
51 using Element = typename GridView::template Codim<0>::Entity;
52
53 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
54 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
55 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
56
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>();
60
61
62public:
63 static constexpr bool isSolDependent = (advectionEnabled && advectionIsSolDependent) ||
64 (diffusionEnabled && diffusionIsSolDependent) ||
65 (heatConductionEnabled && heatConductionIsSolDependent);
66
69 : problemPtr_(&problem) {}
70
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)
90 {
91 // fill the physics-related quantities of the caches
92 if (forceUpdateAll)
93 {
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);
100 }
101 else
102 {
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);
109 }
110 }
111
112private:
113
114 const Problem& problem() const
115 { return *problemPtr_; }
116
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)
124 {
125 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
126 using AdvectionFiller = typename AdvectionType::Cache::Filler;
127
128 // forward to the filler for the advective quantities
129 AdvectionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
130 }
131
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)
139 {
140 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
141 using DiffusionFiller = typename DiffusionType::Cache::Filler;
142 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
143
144 static constexpr int numPhases = ModelTraits::numFluidPhases();
145 static constexpr int numComponents = ModelTraits::numFluidComponents();
146
147 // forward to the filler of the diffusive quantities
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);
152 else
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);
157 }
158
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)
166 {
167 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
168 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
169
170 // forward to the filler of the diffusive quantities
171 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
172 }
173
174 const Problem* problemPtr_;
175};
176
178template<class TypeTag>
179class PorousMediumFluxVariablesCacheFillerImplementation<TypeTag, DiscretizationMethods::CCMpfa>
180{
184 using Element = typename GridView::template Codim<0>::Entity;
186
188 using FVElementGeometry = typename GridGeometry::LocalView;
189 using MpfaHelper = typename GridGeometry::MpfaHelper;
190 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
191 using Extrusion = Extrusion_t<GridGeometry>;
192 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
193 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
194
195 using PrimaryInteractionVolume = GetPropType<TypeTag, Properties::PrimaryInteractionVolume>;
196 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle;
197 using PrimaryLocalFaceData = typename PrimaryInteractionVolume::Traits::LocalFaceData;
198 using SecondaryInteractionVolume = GetPropType<TypeTag, Properties::SecondaryInteractionVolume>;
199 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle;
200 using SecondaryLocalFaceData = typename SecondaryInteractionVolume::Traits::LocalFaceData;
201
202 static constexpr int dim = GridView::dimension;
203 static constexpr int dimWorld = GridView::dimensionworld;
204
205 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
206 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
207 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
208
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>();
212
213public:
217 static constexpr bool isSolDependent = true;
218
221 : problemPtr_(&problem) {}
222
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)
242 {
243 // Set pointers
244 fvGeometryPtr_ = &fvGeometry;
245 elemVolVarsPtr_ = &elemVolVars;
246 const auto& gridGeometry = fvGeometry.gridGeometry();
247
248 // 1. prepare interaction volume (iv)
249 // 2. solve for all transmissibilities and store them in data handles
250 // 3. set pointers to transmissibilities in caches of all the scvfs of the iv
251 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
252 {
253 if (forceUpdateAll)
254 {
255 // create new interaction volume
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);
261
262 // create the corresponding data handle
263 ivDataStorage.secondaryDataHandles.emplace_back();
264 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles.back();
265 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
266
267 // fill the caches for all the scvfs in the interaction volume
268 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
269 }
270 else
271 {
272 // get previously created interaction volume/handle
273 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
274 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes[ivIndexInContainer];
275 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles[ivIndexInContainer];
276 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
277
278 // fill the caches for all the scvfs in the interaction volume
279 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
280 }
281 }
282
283 // primary interaction volume type
284 else
285 {
286 if (forceUpdateAll)
287 {
288 // create new interaction volume
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);
294
295 // create the corresponding data handle
296 ivDataStorage.primaryDataHandles.emplace_back();
297 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles.back();
298 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
299
300 // fill the caches for all the scvfs in the interaction volume
301 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
302 }
303 else
304 {
305 // get previously created interaction volume/handle
306 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
307 primaryIv_ = &ivDataStorage.primaryInteractionVolumes[ivIndexInContainer];
308 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles[ivIndexInContainer];
309 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
310
311 // fill the caches for all the scvfs in the interaction volume
312 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
313 }
314 }
315 }
316
318 const PrimaryInteractionVolume& primaryInteractionVolume() const
319 { return *primaryIv_; }
320
322 const SecondaryInteractionVolume& secondaryInteractionVolume() const
323 { return *secondaryIv_; }
324
326 const PrimaryDataHandle& primaryIvDataHandle() const
327 { return *primaryIvDataHandle_; }
328
330 const SecondaryDataHandle& secondaryIvDataHandle() const
331 { return *secondaryIvDataHandle_; }
332
334 const PrimaryLocalFaceData& primaryIvLocalFaceData() const
335 { return *primaryLocalFaceData_; }
336
338 const SecondaryLocalFaceData& secondaryIvLocalFaceData() const
339 { return *secondaryLocalFaceData_; }
340
341private:
342
343 const Problem& problem() const { return *problemPtr_; }
344 const FVElementGeometry& fvGeometry() const { return *fvGeometryPtr_; }
345 const ElementVolumeVariables& elemVolVars() const { return *elemVolVarsPtr_; }
346
348 template<class FluxVariablesCache, class FluxVarsCacheStorage, class InteractionVolume>
349 void fillCachesInInteractionVolume_(FluxVarsCacheStorage& fluxVarsCacheStorage,
350 InteractionVolume& iv,
351 unsigned int ivIndexInContainer)
352 {
353 // determine if secondary interaction volumes are used here
354 static constexpr bool isSecondary = MpfaHelper::considerSecondaryIVs()
355 && std::is_same_v<InteractionVolume, SecondaryInteractionVolume>;
356
357 // First we update data which are not dependent on the physical processes.
358 // We store pointers to the other flux var caches, so that we have to obtain
359 // this data only once and can use it again in the sub-cache fillers.
360 const auto numGlobalScvfs = iv.localFaceData().size();
361 std::vector<const SubControlVolumeFace*> ivScvfs(numGlobalScvfs);
362 std::vector<FluxVariablesCache*> ivFluxVarCaches(numGlobalScvfs);
363
364 unsigned int i = 0;
365 for (const auto& d : iv.localFaceData())
366 {
367 // obtain the scvf
368 const auto& scvfJ = fvGeometry().scvf(d.gridScvfIndex());
369 ivScvfs[i] = &scvfJ;
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());
375 if (dim < dimWorld)
376 if (d.isOutsideFace())
377 ivFluxVarCaches[i]->setIndexInOutsideFaces(d.scvfLocalOutsideScvfIndex());
378 i++;
379 }
380
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);
387 }
388
390 template<class InteractionVolume, class FluxVariablesCache>
391 void fillAdvection_(InteractionVolume& iv,
392 const std::vector<const SubControlVolumeFace*>& ivScvfs,
393 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
394 {
395 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
396 using AdvectionFiller = typename AdvectionType::Cache::Filler;
397
398 // fill advection caches
399 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
400 {
401 // set pointer to current local face data object
402 // ifs are evaluated at compile time and are optimized away
403 if (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
404 {
405 // we cannot make a distinction, thus we set both pointers
406 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
407 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
408 }
409 else if (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
410 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
411 else
412 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
413
414 // fill this scvfs cache
415 AdvectionFiller::fill(*ivFluxVarCaches[i],
416 problem(),
417 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
418 fvGeometry(),
419 elemVolVars(),
420 *ivScvfs[i],
421 *this);
422 }
423 }
424
426 template<class InteractionVolume, class FluxVariablesCache>
427 void fillDiffusion_(InteractionVolume& iv,
428 const std::vector<const SubControlVolumeFace*>& ivScvfs,
429 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
430 {
431 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
432 using DiffusionFiller = typename DiffusionType::Cache::Filler;
433
434 static constexpr int numPhases = ModelTraits::numFluidPhases();
435 static constexpr int numComponents = ModelTraits::numFluidComponents();
436
437 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
438 {
439 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
440 {
441 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
442 if constexpr (!FluidSystem::isTracerFluidSystem())
443 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
444 continue;
445
446 // fill diffusion caches
447 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
448 {
449 // set pointer to current local face data object
450 // ifs are evaluated at compile time and are optimized away
451 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
452 {
453 // we cannot make a distinction, thus we set both pointers
454 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
455 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
456 }
457 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
458 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
459 else
460 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
461
462 // fill this scvfs cache
463 DiffusionFiller::fill(*ivFluxVarCaches[i],
464 phaseIdx,
465 compIdx,
466 problem(),
467 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
468 fvGeometry(),
469 elemVolVars(),
470 *ivScvfs[i],
471 *this);
472 }
473 }
474 }
475 }
476
478 template<class InteractionVolume, class FluxVariablesCache>
479 void fillHeatConduction_(InteractionVolume& iv,
480 const std::vector<const SubControlVolumeFace*>& ivScvfs,
481 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
482 {
483 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
484 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
485
486 // fill heat conduction caches
487 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
488 {
489 // set pointer to current local face data object
490 // ifs are evaluated at compile time and are optimized away
491 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
492 {
493 // we cannot make a distinction, thus we set both pointers
494 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
495 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
496 }
497 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
498 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
499 else
500 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
501
502 // fill this scvfs cache
503 HeatConductionFiller::fill(*ivFluxVarCaches[i],
504 problem(),
505 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
506 fvGeometry(),
507 elemVolVars(),
508 *ivScvfs[i],
509 *this);
510 }
511 }
512
514 template< class InteractionVolume, class DataHandle>
515 void prepareDataHandle_([[maybe_unused]] InteractionVolume& iv, [[maybe_unused]] DataHandle& handle, [[maybe_unused]] bool forceUpdate)
516 {
517 // (maybe) solve system subject to intrinsic permeability
518 if constexpr (advectionEnabled)
519 {
520 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
521 if constexpr (AdvectionType::discMethod == DiscretizationMethods::ccmpfa)
522 prepareAdvectionHandle_(iv, handle, forceUpdate);
523 }
524
525 // (maybe) solve system subject to diffusion tensors
526 if constexpr (diffusionEnabled)
527 {
528 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
529 if constexpr (DiffusionType::discMethod == DiscretizationMethods::ccmpfa)
530 prepareDiffusionHandles_(iv, handle, forceUpdate);
531 }
532
533 // (maybe) solve system subject to thermal conductivity
534 if constexpr (heatConductionEnabled)
535 {
536 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
537 if constexpr (HeatConductionType::discMethod == DiscretizationMethods::ccmpfa)
538 prepareHeatConductionHandle_(iv, handle, forceUpdate);
539 }
540 }
541
543 template<class InteractionVolume, class DataHandle>
544 void prepareAdvectionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
545 {
546 // get instance of the interaction volume-local assembler
547 using Traits = typename InteractionVolume::Traits;
548 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
549 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
550
551 // lambda to obtain the permeability tensor
552 auto getK = [] (const auto& volVars) { return volVars.permeability(); };
553
554 // Assemble T only if permeability is sol-dependent or if update is forced
555 if (forceUpdateAll || advectionIsSolDependent)
556 localAssembler.assembleMatrices(handle.advectionHandle(), iv, getK);
557
558 // assemble pressure vectors
559 for (unsigned int pIdx = 0; pIdx < ModelTraits::numFluidPhases(); ++pIdx)
560 {
561 // set context in handle
562 handle.advectionHandle().setPhaseIndex(pIdx);
563
564 // maybe (re-)assemble gravity contribution vector
565 auto getRho = [pIdx] (const auto& volVars) { return volVars.density(pIdx); };
566 static const bool enableGravity = getParamFromGroup<bool>(problem().paramGroup(), "Problem.EnableGravity");
567 if (enableGravity)
568 localAssembler.assembleGravity(handle.advectionHandle(), iv, getRho);
569
570 // reassemble pressure vector
571 auto getPressure = [pIdx] (const auto& volVars) { return volVars.pressure(pIdx); };
572 localAssembler.assembleU(handle.advectionHandle(), iv, getPressure);
573 }
574 }
575
577 template<class InteractionVolume, class DataHandle>
578 void prepareDiffusionHandles_(InteractionVolume& iv,
579 DataHandle& handle,
580 bool forceUpdateAll)
581 {
582 for (unsigned int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
583 {
584 for (unsigned int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
585 {
586 // skip main component
587 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
588 if constexpr (!FluidSystem::isTracerFluidSystem())
589 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
590 continue;
591
592 // fill data in the handle
593 handle.diffusionHandle().setPhaseIndex(phaseIdx);
594 handle.diffusionHandle().setComponentIndex(compIdx);
595
596 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
597
598 // get instance of the interaction volume-local assembler
599 using Traits = typename InteractionVolume::Traits;
600 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
601 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
602
603 // maybe (re-)assemble matrices
604 if (forceUpdateAll || diffusionIsSolDependent)
605 {
606 // lambda to obtain diffusion coefficient
607 const auto getD = [&](const auto& volVars)
608 {
609 if constexpr (FluidSystem::isTracerFluidSystem())
610 return volVars.effectiveDiffusionCoefficient(0, 0, compIdx);
611 else
612 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
613 };
614
615 // Effective diffusion coefficients might be zero if saturation = 0.
616 // Compute epsilon to detect obsolete rows in the iv-local matrices during assembly
617 static const auto zeroD = getParamFromGroup<Scalar>(
618 problem().paramGroup(),
619 "Mpfa.ZeroEffectiveDiffusionCoefficientThreshold",
620 1e-16
621 );
622
623 // compute a representative transmissibility for this interaction volume, using
624 // the threshold for zero diffusion coefficients, and use this as epsilon
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];
628 const auto eps = Extrusion::area(fvGeometry(), scvf)*computeTpfaTransmissibility(
629 fvGeometry(), scvf, scv, zeroD, vv.extrusionFactor()
630 );
631
632 localAssembler.assembleMatrices(handle.diffusionHandle(), iv, getD, eps);
633 }
634
635 // assemble vector of mole fractions
636 auto getMassOrMoleFraction = [phaseIdx, compIdx] (const auto& volVars)
637 {
638 return (DiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged) ? volVars.massFraction(phaseIdx, compIdx) :
639 volVars.moleFraction(phaseIdx, compIdx);
640 };
641
642 localAssembler.assembleU(handle.diffusionHandle(), iv, getMassOrMoleFraction);
643 }
644 }
645 }
646
648 template<class InteractionVolume, class DataHandle>
649 void prepareHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
650 {
651 // get instance of the interaction volume-local assembler
652 using Traits = typename InteractionVolume::Traits;
653 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
654 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
655
656 // lambda to obtain the effective thermal conductivity
657 auto getLambda = [] (const auto& volVars) { return volVars.effectiveThermalConductivity(); };
658
659 // maybe (re-)assemble matrices
660 if (forceUpdateAll || heatConductionIsSolDependent)
661 localAssembler.assembleMatrices(handle.heatConductionHandle(), iv, getLambda);
662
663 // assemble vector of temperatures
664 auto getTemperature = [] (const auto& volVars) { return volVars.temperature(); };
665 localAssembler.assembleU(handle.heatConductionHandle(), iv, getTemperature);
666 }
667
668 const Problem* problemPtr_;
669 const FVElementGeometry* fvGeometryPtr_;
670 const ElementVolumeVariables* elemVolVarsPtr_;
671
672 // We store pointers to an inner and a boundary interaction volume.
673 // These are updated during the filling of the caches and the
674 // physics-related caches have access to them
675 PrimaryInteractionVolume* primaryIv_;
676 SecondaryInteractionVolume* secondaryIv_;
677
678 // pointer to the current interaction volume data handle
679 PrimaryDataHandle* primaryIvDataHandle_;
680 SecondaryDataHandle* secondaryIvDataHandle_;
681
682 // We do an interaction volume-wise filling of the caches
683 // While filling, we store a pointer to the current localScvf
684 // face data object of the IV so that the individual caches
685 // can access it and don't have to retrieve it again
686 const PrimaryLocalFaceData* primaryLocalFaceData_;
687 const SecondaryLocalFaceData* secondaryLocalFaceData_;
688};
689
690} // end namespace Dumux
691
692#endif
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
Definition: adapt.hh:17
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.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...