3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
25#define DUMUX_POROUSMEDIUM_FLUXVARIABLESCACHE_FILLER_HH
26
29
34
35namespace Dumux {
36
37// forward declaration
38template<class TypeTag, DiscretizationMethod discMethod>
40
47template<class TypeTag>
49
51template<class TypeTag>
53{
57 using GridView = typename GridGeometry::GridView;
58 using FVElementGeometry = typename GridGeometry::LocalView;
59 using SubControlVolume = typename GridGeometry::SubControlVolume;
60 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
61 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
62
63 using Element = typename GridView::template Codim<0>::Entity;
64
65 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
66 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
67 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
68
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>();
72
73
74public:
75 static constexpr bool isSolDependent = (advectionEnabled && advectionIsSolDependent) ||
76 (diffusionEnabled && diffusionIsSolDependent) ||
77 (heatConductionEnabled && heatConductionIsSolDependent);
78
81 : problemPtr_(&problem) {}
82
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)
102 {
103 // fill the physics-related quantities of the caches
104 if (forceUpdateAll)
105 {
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);
112 }
113 else
114 {
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);
121 }
122 }
123
124private:
125
126 const Problem& problem() const
127 { return *problemPtr_; }
128
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)
136 {
137 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
138 using AdvectionFiller = typename AdvectionType::Cache::Filler;
139
140 // forward to the filler for the advective quantities
141 AdvectionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
142 }
143
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)
151 {
152 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
153 using DiffusionFiller = typename DiffusionType::Cache::Filler;
154 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
155
156 static constexpr int numPhases = ModelTraits::numFluidPhases();
157 static constexpr int numComponents = ModelTraits::numFluidComponents();
158
159 // forward to the filler of the diffusive quantities
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);
164 else
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);
169 }
170
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)
178 {
179 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
180 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
181
182 // forward to the filler of the diffusive quantities
183 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
184 }
185
186 const Problem* problemPtr_;
187};
188
190template<class TypeTag>
192{
196 using Element = typename GridView::template Codim<0>::Entity;
197
199 using FVElementGeometry = typename GridGeometry::LocalView;
200 using MpfaHelper = typename GridGeometry::MpfaHelper;
201 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
202 using Extrusion = Extrusion_t<GridGeometry>;
203 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
204 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
205
206 using PrimaryInteractionVolume = GetPropType<TypeTag, Properties::PrimaryInteractionVolume>;
207 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle;
208 using PrimaryLocalFaceData = typename PrimaryInteractionVolume::Traits::LocalFaceData;
209 using SecondaryInteractionVolume = GetPropType<TypeTag, Properties::SecondaryInteractionVolume>;
210 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle;
211 using SecondaryLocalFaceData = typename SecondaryInteractionVolume::Traits::LocalFaceData;
212
213 static constexpr int dim = GridView::dimension;
214 static constexpr int dimWorld = GridView::dimensionworld;
215
216 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
217 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
218 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
219
220 static constexpr bool advectionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
221 static constexpr bool diffusionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
222 static constexpr bool heatConductionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
223
224public:
228 static constexpr bool isSolDependent = true;
229
232 : problemPtr_(&problem) {}
233
246 template<class FluxVarsCacheStorage, class FluxVariablesCache, class IVDataStorage>
247 void fill(FluxVarsCacheStorage& fluxVarsCacheStorage,
248 FluxVariablesCache& scvfFluxVarsCache,
249 IVDataStorage& ivDataStorage,
250 const Element& element,
251 const FVElementGeometry& fvGeometry,
252 const ElementVolumeVariables& elemVolVars,
253 const SubControlVolumeFace& scvf,
254 bool forceUpdateAll = false)
255 {
256 // Set pointers
257 elementPtr_ = &element;
258 fvGeometryPtr_ = &fvGeometry;
259 elemVolVarsPtr_ = &elemVolVars;
260 const auto& gridGeometry = fvGeometry.gridGeometry();
261
262 // 1. prepare interaction volume (iv)
263 // 2. solve for all transmissibilities and store them in data handles
264 // 3. set pointers to transmissibilities in caches of all the scvfs of the iv
265 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
266 {
267 if (forceUpdateAll)
268 {
269 // create new interaction volume
270 const auto ivIndexInContainer = ivDataStorage.secondaryInteractionVolumes.size();
271 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().secondaryIndexSet(scvf);
272 ivDataStorage.secondaryInteractionVolumes.emplace_back();
273 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes.back();
274 secondaryIv_->bind(indexSet, problem(), fvGeometry);
275
276 // create the corresponding data handle
277 ivDataStorage.secondaryDataHandles.emplace_back();
278 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles.back();
279 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
280
281 // fill the caches for all the scvfs in the interaction volume
282 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
283 }
284 else
285 {
286 // get previously created interaction volume/handle
287 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
288 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes[ivIndexInContainer];
289 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles[ivIndexInContainer];
290 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
291
292 // fill the caches for all the scvfs in the interaction volume
293 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
294 }
295 }
296
297 // primary interaction volume type
298 else
299 {
300 if (forceUpdateAll)
301 {
302 // create new interaction volume
303 const auto ivIndexInContainer = ivDataStorage.primaryInteractionVolumes.size();
304 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().primaryIndexSet(scvf);
305 ivDataStorage.primaryInteractionVolumes.emplace_back();
306 primaryIv_ = &ivDataStorage.primaryInteractionVolumes.back();
307 primaryIv_->bind(indexSet, problem(), fvGeometry);
308
309 // create the corresponding data handle
310 ivDataStorage.primaryDataHandles.emplace_back();
311 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles.back();
312 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
313
314 // fill the caches for all the scvfs in the interaction volume
315 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
316 }
317 else
318 {
319 // get previously created interaction volume/handle
320 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
321 primaryIv_ = &ivDataStorage.primaryInteractionVolumes[ivIndexInContainer];
322 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles[ivIndexInContainer];
323 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
324
325 // fill the caches for all the scvfs in the interaction volume
326 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
327 }
328 }
329 }
330
332 const PrimaryInteractionVolume& primaryInteractionVolume() const
333 { return *primaryIv_; }
334
336 const SecondaryInteractionVolume& secondaryInteractionVolume() const
337 { return *secondaryIv_; }
338
340 const PrimaryDataHandle& primaryIvDataHandle() const
341 { return *primaryIvDataHandle_; }
342
344 const SecondaryDataHandle& secondaryIvDataHandle() const
345 { return *secondaryIvDataHandle_; }
346
348 const PrimaryLocalFaceData& primaryIvLocalFaceData() const
349 { return *primaryLocalFaceData_; }
350
352 const SecondaryLocalFaceData& secondaryIvLocalFaceData() const
353 { return *secondaryLocalFaceData_; }
354
355private:
356
357 const Problem& problem() const { return *problemPtr_; }
358 const Element& element() const { return *elementPtr_; }
359 const FVElementGeometry& fvGeometry() const { return *fvGeometryPtr_; }
360 const ElementVolumeVariables& elemVolVars() const { return *elemVolVarsPtr_; }
361
363 template<class FluxVariablesCache, class FluxVarsCacheStorage, class InteractionVolume>
364 void fillCachesInInteractionVolume_(FluxVarsCacheStorage& fluxVarsCacheStorage,
365 InteractionVolume& iv,
366 unsigned int ivIndexInContainer)
367 {
368 // determine if secondary interaction volumes are used here
369 static constexpr bool isSecondary = MpfaHelper::considerSecondaryIVs()
370 && std::is_same_v<InteractionVolume, SecondaryInteractionVolume>;
371
372 // First we upate data which are not dependent on the physical processes.
373 // We store pointers to the other flux var caches, so that we have to obtain
374 // this data only once and can use it again in the sub-cache fillers.
375 const auto numGlobalScvfs = iv.localFaceData().size();
376 std::vector<const SubControlVolumeFace*> ivScvfs(numGlobalScvfs);
377 std::vector<FluxVariablesCache*> ivFluxVarCaches(numGlobalScvfs);
378
379 unsigned int i = 0;
380 for (const auto& d : iv.localFaceData())
381 {
382 // obtain the scvf
383 const auto& scvfJ = fvGeometry().scvf(d.gridScvfIndex());
384 ivScvfs[i] = &scvfJ;
385 ivFluxVarCaches[i] = &fluxVarsCacheStorage[scvfJ];
386 ivFluxVarCaches[i]->setIvIndexInContainer(ivIndexInContainer);
387 ivFluxVarCaches[i]->setUpdateStatus(true);
388 ivFluxVarCaches[i]->setSecondaryIvUsage(isSecondary);
389 ivFluxVarCaches[i]->setIvLocalFaceIndex(d.ivLocalScvfIndex());
390 if (dim < dimWorld)
391 if (d.isOutsideFace())
392 ivFluxVarCaches[i]->setIndexInOutsideFaces(d.scvfLocalOutsideScvfIndex());
393 i++;
394 }
395
396 if constexpr (advectionEnabled)
397 fillAdvection_(iv, ivScvfs, ivFluxVarCaches);
398 if constexpr (diffusionEnabled)
399 fillDiffusion_(iv, ivScvfs, ivFluxVarCaches);
400 if constexpr (heatConductionEnabled)
401 fillHeatConduction_(iv, ivScvfs, ivFluxVarCaches);
402 }
403
405 template<class InteractionVolume, class FluxVariablesCache>
406 void fillAdvection_(InteractionVolume& iv,
407 const std::vector<const SubControlVolumeFace*>& ivScvfs,
408 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
409 {
410 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
411 using AdvectionFiller = typename AdvectionType::Cache::Filler;
412
413 // fill advection caches
414 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
415 {
416 // set pointer to current local face data object
417 // ifs are evaluated at compile time and are optimized away
418 if (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
419 {
420 // we cannot make a disctinction, thus we set both pointers
421 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
422 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
423 }
424 else if (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
425 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
426 else
427 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
428
429 // fill this scvfs cache
430 AdvectionFiller::fill(*ivFluxVarCaches[i],
431 problem(),
432 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
433 fvGeometry(),
434 elemVolVars(),
435 *ivScvfs[i],
436 *this);
437 }
438 }
439
441 template<class InteractionVolume, class FluxVariablesCache>
442 void fillDiffusion_(InteractionVolume& iv,
443 const std::vector<const SubControlVolumeFace*>& ivScvfs,
444 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
445 {
446 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
447 using DiffusionFiller = typename DiffusionType::Cache::Filler;
448
449 static constexpr int numPhases = ModelTraits::numFluidPhases();
450 static constexpr int numComponents = ModelTraits::numFluidComponents();
451
452 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
453 {
454 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
455 {
456 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
457 if constexpr (!FluidSystem::isTracerFluidSystem())
458 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
459 continue;
460
461 // fill diffusion caches
462 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
463 {
464 // set pointer to current local face data object
465 // ifs are evaluated at compile time and are optimized away
466 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
467 {
468 // we cannot make a disctinction, thus we set both pointers
469 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
470 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
471 }
472 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
473 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
474 else
475 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
476
477 // fill this scvfs cache
478 DiffusionFiller::fill(*ivFluxVarCaches[i],
479 phaseIdx,
480 compIdx,
481 problem(),
482 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
483 fvGeometry(),
484 elemVolVars(),
485 *ivScvfs[i],
486 *this);
487 }
488 }
489 }
490 }
491
493 template<class InteractionVolume, class FluxVariablesCache>
494 void fillHeatConduction_(InteractionVolume& iv,
495 const std::vector<const SubControlVolumeFace*>& ivScvfs,
496 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
497 {
498 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
499 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
500
501 // fill heat conduction caches
502 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
503 {
504 // set pointer to current local face data object
505 // ifs are evaluated at compile time and are optimized away
506 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
507 {
508 // we cannot make a disctinction, thus we set both pointers
509 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
510 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
511 }
512 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
513 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
514 else
515 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
516
517 // fill this scvfs cache
518 HeatConductionFiller::fill(*ivFluxVarCaches[i],
519 problem(),
520 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
521 fvGeometry(),
522 elemVolVars(),
523 *ivScvfs[i],
524 *this);
525 }
526 }
527
529 template< class InteractionVolume, class DataHandle>
530 void prepareDataHandle_([[maybe_unused]] InteractionVolume& iv, [[maybe_unused]] DataHandle& handle, [[maybe_unused]] bool forceUpdate)
531 {
532 // (maybe) solve system subject to intrinsic permeability
533 if constexpr (advectionEnabled)
534 {
535 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
536 if constexpr (AdvectionType::discMethod == DiscretizationMethod::ccmpfa)
537 prepareAdvectionHandle_(iv, handle, forceUpdate);
538 }
539
540 // (maybe) solve system subject to diffusion tensors
541 if constexpr (diffusionEnabled)
542 {
543 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
544 if constexpr (DiffusionType::discMethod == DiscretizationMethod::ccmpfa)
545 prepareDiffusionHandles_(iv, handle, forceUpdate);
546 }
547
548 // (maybe) solve system subject to thermal conductivity
549 if constexpr (heatConductionEnabled)
550 {
551 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
552 if constexpr (HeatConductionType::discMethod == DiscretizationMethod::ccmpfa)
553 prepareHeatConductionHandle_(iv, handle, forceUpdate);
554 }
555 }
556
558 template<class InteractionVolume, class DataHandle>
559 void prepareAdvectionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
560 {
561 // get instance of the interaction volume-local assembler
562 using Traits = typename InteractionVolume::Traits;
563 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
564 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
565
566 // lambda to obtain the permeability tensor
567 auto getK = [] (const auto& volVars) { return volVars.permeability(); };
568
569 // Assemble T only if permeability is sol-dependent or if update is forced
570 if (forceUpdateAll || advectionIsSolDependent)
571 localAssembler.assembleMatrices(handle.advectionHandle(), iv, getK);
572
573 // assemble pressure vectors
574 for (unsigned int pIdx = 0; pIdx < ModelTraits::numFluidPhases(); ++pIdx)
575 {
576 // set context in handle
577 handle.advectionHandle().setPhaseIndex(pIdx);
578
579 // maybe (re-)assemble gravity contribution vector
580 auto getRho = [pIdx] (const auto& volVars) { return volVars.density(pIdx); };
581 static const bool enableGravity = getParamFromGroup<bool>(problem().paramGroup(), "Problem.EnableGravity");
582 if (enableGravity)
583 localAssembler.assembleGravity(handle.advectionHandle(), iv, getRho);
584
585 // reassemble pressure vector
586 auto getPressure = [pIdx] (const auto& volVars) { return volVars.pressure(pIdx); };
587 localAssembler.assembleU(handle.advectionHandle(), iv, getPressure);
588 }
589 }
590
592 template<class InteractionVolume, class DataHandle>
593 void prepareDiffusionHandles_(InteractionVolume& iv,
594 DataHandle& handle,
595 bool forceUpdateAll)
596 {
597 for (unsigned int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
598 {
599 for (unsigned int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
600 {
601 // skip main component
602 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
603 if constexpr (!FluidSystem::isTracerFluidSystem())
604 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
605 continue;
606
607 // fill data in the handle
608 handle.diffusionHandle().setPhaseIndex(phaseIdx);
609 handle.diffusionHandle().setComponentIndex(compIdx);
610
611 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
612
613 // get instance of the interaction volume-local assembler
614 using Traits = typename InteractionVolume::Traits;
615 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
616 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
617
618 // maybe (re-)assemble matrices
619 if (forceUpdateAll || diffusionIsSolDependent)
620 {
621 // lambda to obtain diffusion coefficient
622 const auto getD = [&](const auto& volVars)
623 {
624 if constexpr (FluidSystem::isTracerFluidSystem())
625 return volVars.effectiveDiffusionCoefficient(0, 0, compIdx);
626 else
627 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
628 };
629
630 // Effective diffusion coefficients might get zero if saturation = 0.
631 // Compute epsilon to detect obsolete rows in the iv-local matrices during assembly
632 const auto& scv = *scvs(fvGeometry()).begin();
633 const auto& scvf = *scvfs(fvGeometry()).begin();
634 const auto& vv = elemVolVars()[scv];
635 const auto D = [&] ()
636 {
637 // diffusion coefficients below 1e-20 are treated as zeroes!!
638 using std::max;
639 if constexpr (!FluidSystem::isTracerFluidSystem())
640 return max(1e-20, vv.diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx));
641 else
642 return max(1e-20, vv.diffusionCoefficient(0, 0, compIdx));
643 } ();
644
645 auto eps = 1e-7*computeTpfaTransmissibility(scvf, scv, D, vv.extrusionFactor())*Extrusion::area(scvf);
646 localAssembler.assembleMatrices(handle.diffusionHandle(), iv, getD, eps);
647 }
648
649 // assemble vector of mole fractions
650 auto getMassOrMoleFraction = [phaseIdx, compIdx] (const auto& volVars)
651 {
652 return (DiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged) ? volVars.massFraction(phaseIdx, compIdx) :
653 volVars.moleFraction(phaseIdx, compIdx);
654 };
655
656 localAssembler.assembleU(handle.diffusionHandle(), iv, getMassOrMoleFraction);
657 }
658 }
659 }
660
662 template<class InteractionVolume, class DataHandle>
663 void prepareHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
664 {
665 // get instance of the interaction volume-local assembler
666 using Traits = typename InteractionVolume::Traits;
667 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
668 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
669
670 // lambda to obtain the effective thermal conductivity
671 auto getLambda = [] (const auto& volVars) { return volVars.effectiveThermalConductivity(); };
672
673 // maybe (re-)assemble matrices
674 if (forceUpdateAll || heatConductionIsSolDependent)
675 localAssembler.assembleMatrices(handle.heatConductionHandle(), iv, getLambda);
676
677 // assemble vector of temperatures
678 auto getTemperature = [] (const auto& volVars) { return volVars.temperature(); };
679 localAssembler.assembleU(handle.heatConductionHandle(), iv, getTemperature);
680 }
681
682 const Problem* problemPtr_;
683 const Element* elementPtr_;
684 const FVElementGeometry* fvGeometryPtr_;
685 const ElementVolumeVariables* elemVolVarsPtr_;
686
687 // We store pointers to an inner and a boundary interaction volume.
688 // These are updated during the filling of the caches and the
689 // physics-related caches have access to them
690 PrimaryInteractionVolume* primaryIv_;
691 SecondaryInteractionVolume* secondaryIv_;
692
693 // pointer to the current interaction volume data handle
694 PrimaryDataHandle* primaryIvDataHandle_;
695 SecondaryDataHandle* secondaryIvDataHandle_;
696
697 // We do an interaction volume-wise filling of the caches
698 // While filling, we store a pointer to the current localScvf
699 // face data object of the IV so that the individual caches
700 // can access it and don't have to retrieve it again
701 const PrimaryLocalFaceData* primaryLocalFaceData_;
702 const SecondaryLocalFaceData* secondaryLocalFaceData_;
703};
704
705} // end namespace Dumux
706
707#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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
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
Definition: adapt.hh:29
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 (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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 PrimaryLocalFaceData & primaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:348
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: fluxvariablescachefiller.hh:247
const PrimaryInteractionVolume & primaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:332
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets problem pointer.
Definition: fluxvariablescachefiller.hh:231
const PrimaryDataHandle & primaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:340
const SecondaryDataHandle & secondaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:344
const SecondaryLocalFaceData & secondaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:352
const SecondaryInteractionVolume & secondaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:336
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...