3.5-git
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, class DiscretizationMethod>
40
47template<class TypeTag>
49
51template<class TypeTag>
52class PorousMediumFluxVariablesCacheFillerImplementation<TypeTag, DiscretizationMethods::CCTpfa>
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>
191class PorousMediumFluxVariablesCacheFillerImplementation<TypeTag, DiscretizationMethods::CCMpfa>
192{
196 using Element = typename GridView::template Codim<0>::Entity;
198
200 using FVElementGeometry = typename GridGeometry::LocalView;
201 using MpfaHelper = typename GridGeometry::MpfaHelper;
202 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
203 using Extrusion = Extrusion_t<GridGeometry>;
204 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
205 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
206
207 using PrimaryInteractionVolume = GetPropType<TypeTag, Properties::PrimaryInteractionVolume>;
208 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle;
209 using PrimaryLocalFaceData = typename PrimaryInteractionVolume::Traits::LocalFaceData;
210 using SecondaryInteractionVolume = GetPropType<TypeTag, Properties::SecondaryInteractionVolume>;
211 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle;
212 using SecondaryLocalFaceData = typename SecondaryInteractionVolume::Traits::LocalFaceData;
213
214 static constexpr int dim = GridView::dimension;
215 static constexpr int dimWorld = GridView::dimensionworld;
216
217 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
218 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
219 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
220
221 static constexpr bool advectionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
222 static constexpr bool diffusionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
223 static constexpr bool heatConductionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
224
225public:
229 static constexpr bool isSolDependent = true;
230
233 : problemPtr_(&problem) {}
234
246 template<class FluxVarsCacheStorage, class FluxVariablesCache, class IVDataStorage>
247 void fill(FluxVarsCacheStorage& fluxVarsCacheStorage,
248 FluxVariablesCache& scvfFluxVarsCache,
249 IVDataStorage& ivDataStorage,
250 const FVElementGeometry& fvGeometry,
251 const ElementVolumeVariables& elemVolVars,
252 const SubControlVolumeFace& scvf,
253 bool forceUpdateAll = false)
254 {
255 // Set pointers
256 fvGeometryPtr_ = &fvGeometry;
257 elemVolVarsPtr_ = &elemVolVars;
258 const auto& gridGeometry = fvGeometry.gridGeometry();
259
260 // 1. prepare interaction volume (iv)
261 // 2. solve for all transmissibilities and store them in data handles
262 // 3. set pointers to transmissibilities in caches of all the scvfs of the iv
263 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
264 {
265 if (forceUpdateAll)
266 {
267 // create new interaction volume
268 const auto ivIndexInContainer = ivDataStorage.secondaryInteractionVolumes.size();
269 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().secondaryIndexSet(scvf);
270 ivDataStorage.secondaryInteractionVolumes.emplace_back();
271 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes.back();
272 secondaryIv_->bind(indexSet, problem(), fvGeometry);
273
274 // create the corresponding data handle
275 ivDataStorage.secondaryDataHandles.emplace_back();
276 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles.back();
277 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
278
279 // fill the caches for all the scvfs in the interaction volume
280 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
281 }
282 else
283 {
284 // get previously created interaction volume/handle
285 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
286 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes[ivIndexInContainer];
287 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles[ivIndexInContainer];
288 prepareDataHandle_(*secondaryIv_, *secondaryIvDataHandle_, forceUpdateAll);
289
290 // fill the caches for all the scvfs in the interaction volume
291 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, ivIndexInContainer);
292 }
293 }
294
295 // primary interaction volume type
296 else
297 {
298 if (forceUpdateAll)
299 {
300 // create new interaction volume
301 const auto ivIndexInContainer = ivDataStorage.primaryInteractionVolumes.size();
302 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().primaryIndexSet(scvf);
303 ivDataStorage.primaryInteractionVolumes.emplace_back();
304 primaryIv_ = &ivDataStorage.primaryInteractionVolumes.back();
305 primaryIv_->bind(indexSet, problem(), fvGeometry);
306
307 // create the corresponding data handle
308 ivDataStorage.primaryDataHandles.emplace_back();
309 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles.back();
310 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
311
312 // fill the caches for all the scvfs in the interaction volume
313 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
314 }
315 else
316 {
317 // get previously created interaction volume/handle
318 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
319 primaryIv_ = &ivDataStorage.primaryInteractionVolumes[ivIndexInContainer];
320 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles[ivIndexInContainer];
321 prepareDataHandle_(*primaryIv_, *primaryIvDataHandle_, forceUpdateAll);
322
323 // fill the caches for all the scvfs in the interaction volume
324 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, ivIndexInContainer);
325 }
326 }
327 }
328
330 const PrimaryInteractionVolume& primaryInteractionVolume() const
331 { return *primaryIv_; }
332
334 const SecondaryInteractionVolume& secondaryInteractionVolume() const
335 { return *secondaryIv_; }
336
338 const PrimaryDataHandle& primaryIvDataHandle() const
339 { return *primaryIvDataHandle_; }
340
342 const SecondaryDataHandle& secondaryIvDataHandle() const
343 { return *secondaryIvDataHandle_; }
344
346 const PrimaryLocalFaceData& primaryIvLocalFaceData() const
347 { return *primaryLocalFaceData_; }
348
350 const SecondaryLocalFaceData& secondaryIvLocalFaceData() const
351 { return *secondaryLocalFaceData_; }
352
353private:
354
355 const Problem& problem() const { return *problemPtr_; }
356 const FVElementGeometry& fvGeometry() const { return *fvGeometryPtr_; }
357 const ElementVolumeVariables& elemVolVars() const { return *elemVolVarsPtr_; }
358
360 template<class FluxVariablesCache, class FluxVarsCacheStorage, class InteractionVolume>
361 void fillCachesInInteractionVolume_(FluxVarsCacheStorage& fluxVarsCacheStorage,
362 InteractionVolume& iv,
363 unsigned int ivIndexInContainer)
364 {
365 // determine if secondary interaction volumes are used here
366 static constexpr bool isSecondary = MpfaHelper::considerSecondaryIVs()
367 && std::is_same_v<InteractionVolume, SecondaryInteractionVolume>;
368
369 // First we upate data which are not dependent on the physical processes.
370 // We store pointers to the other flux var caches, so that we have to obtain
371 // this data only once and can use it again in the sub-cache fillers.
372 const auto numGlobalScvfs = iv.localFaceData().size();
373 std::vector<const SubControlVolumeFace*> ivScvfs(numGlobalScvfs);
374 std::vector<FluxVariablesCache*> ivFluxVarCaches(numGlobalScvfs);
375
376 unsigned int i = 0;
377 for (const auto& d : iv.localFaceData())
378 {
379 // obtain the scvf
380 const auto& scvfJ = fvGeometry().scvf(d.gridScvfIndex());
381 ivScvfs[i] = &scvfJ;
382 ivFluxVarCaches[i] = &fluxVarsCacheStorage[scvfJ];
383 ivFluxVarCaches[i]->setIvIndexInContainer(ivIndexInContainer);
384 ivFluxVarCaches[i]->setUpdateStatus(true);
385 ivFluxVarCaches[i]->setSecondaryIvUsage(isSecondary);
386 ivFluxVarCaches[i]->setIvLocalFaceIndex(d.ivLocalScvfIndex());
387 if (dim < dimWorld)
388 if (d.isOutsideFace())
389 ivFluxVarCaches[i]->setIndexInOutsideFaces(d.scvfLocalOutsideScvfIndex());
390 i++;
391 }
392
393 if constexpr (advectionEnabled)
394 fillAdvection_(iv, ivScvfs, ivFluxVarCaches);
395 if constexpr (diffusionEnabled)
396 fillDiffusion_(iv, ivScvfs, ivFluxVarCaches);
397 if constexpr (heatConductionEnabled)
398 fillHeatConduction_(iv, ivScvfs, ivFluxVarCaches);
399 }
400
402 template<class InteractionVolume, class FluxVariablesCache>
403 void fillAdvection_(InteractionVolume& iv,
404 const std::vector<const SubControlVolumeFace*>& ivScvfs,
405 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
406 {
407 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
408 using AdvectionFiller = typename AdvectionType::Cache::Filler;
409
410 // fill advection caches
411 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
412 {
413 // set pointer to current local face data object
414 // ifs are evaluated at compile time and are optimized away
415 if (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
416 {
417 // we cannot make a disctinction, thus we set both pointers
418 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
419 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
420 }
421 else if (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
422 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
423 else
424 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
425
426 // fill this scvfs cache
427 AdvectionFiller::fill(*ivFluxVarCaches[i],
428 problem(),
429 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
430 fvGeometry(),
431 elemVolVars(),
432 *ivScvfs[i],
433 *this);
434 }
435 }
436
438 template<class InteractionVolume, class FluxVariablesCache>
439 void fillDiffusion_(InteractionVolume& iv,
440 const std::vector<const SubControlVolumeFace*>& ivScvfs,
441 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
442 {
443 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
444 using DiffusionFiller = typename DiffusionType::Cache::Filler;
445
446 static constexpr int numPhases = ModelTraits::numFluidPhases();
447 static constexpr int numComponents = ModelTraits::numFluidComponents();
448
449 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
450 {
451 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
452 {
453 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
454 if constexpr (!FluidSystem::isTracerFluidSystem())
455 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
456 continue;
457
458 // fill diffusion caches
459 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
460 {
461 // set pointer to current local face data object
462 // ifs are evaluated at compile time and are optimized away
463 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
464 {
465 // we cannot make a disctinction, thus we set both pointers
466 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
467 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
468 }
469 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
470 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
471 else
472 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
473
474 // fill this scvfs cache
475 DiffusionFiller::fill(*ivFluxVarCaches[i],
476 phaseIdx,
477 compIdx,
478 problem(),
479 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
480 fvGeometry(),
481 elemVolVars(),
482 *ivScvfs[i],
483 *this);
484 }
485 }
486 }
487 }
488
490 template<class InteractionVolume, class FluxVariablesCache>
491 void fillHeatConduction_(InteractionVolume& iv,
492 const std::vector<const SubControlVolumeFace*>& ivScvfs,
493 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
494 {
495 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
496 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
497
498 // fill heat conduction caches
499 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
500 {
501 // set pointer to current local face data object
502 // ifs are evaluated at compile time and are optimized away
503 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
504 {
505 // we cannot make a disctinction, thus we set both pointers
506 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
507 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
508 }
509 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
510 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
511 else
512 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
513
514 // fill this scvfs cache
515 HeatConductionFiller::fill(*ivFluxVarCaches[i],
516 problem(),
517 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
518 fvGeometry(),
519 elemVolVars(),
520 *ivScvfs[i],
521 *this);
522 }
523 }
524
526 template< class InteractionVolume, class DataHandle>
527 void prepareDataHandle_([[maybe_unused]] InteractionVolume& iv, [[maybe_unused]] DataHandle& handle, [[maybe_unused]] bool forceUpdate)
528 {
529 // (maybe) solve system subject to intrinsic permeability
530 if constexpr (advectionEnabled)
531 {
532 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
533 if constexpr (AdvectionType::discMethod == DiscretizationMethods::ccmpfa)
534 prepareAdvectionHandle_(iv, handle, forceUpdate);
535 }
536
537 // (maybe) solve system subject to diffusion tensors
538 if constexpr (diffusionEnabled)
539 {
540 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
541 if constexpr (DiffusionType::discMethod == DiscretizationMethods::ccmpfa)
542 prepareDiffusionHandles_(iv, handle, forceUpdate);
543 }
544
545 // (maybe) solve system subject to thermal conductivity
546 if constexpr (heatConductionEnabled)
547 {
548 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
549 if constexpr (HeatConductionType::discMethod == DiscretizationMethods::ccmpfa)
550 prepareHeatConductionHandle_(iv, handle, forceUpdate);
551 }
552 }
553
555 template<class InteractionVolume, class DataHandle>
556 void prepareAdvectionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
557 {
558 // get instance of the interaction volume-local assembler
559 using Traits = typename InteractionVolume::Traits;
560 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
561 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
562
563 // lambda to obtain the permeability tensor
564 auto getK = [] (const auto& volVars) { return volVars.permeability(); };
565
566 // Assemble T only if permeability is sol-dependent or if update is forced
567 if (forceUpdateAll || advectionIsSolDependent)
568 localAssembler.assembleMatrices(handle.advectionHandle(), iv, getK);
569
570 // assemble pressure vectors
571 for (unsigned int pIdx = 0; pIdx < ModelTraits::numFluidPhases(); ++pIdx)
572 {
573 // set context in handle
574 handle.advectionHandle().setPhaseIndex(pIdx);
575
576 // maybe (re-)assemble gravity contribution vector
577 auto getRho = [pIdx] (const auto& volVars) { return volVars.density(pIdx); };
578 static const bool enableGravity = getParamFromGroup<bool>(problem().paramGroup(), "Problem.EnableGravity");
579 if (enableGravity)
580 localAssembler.assembleGravity(handle.advectionHandle(), iv, getRho);
581
582 // reassemble pressure vector
583 auto getPressure = [pIdx] (const auto& volVars) { return volVars.pressure(pIdx); };
584 localAssembler.assembleU(handle.advectionHandle(), iv, getPressure);
585 }
586 }
587
589 template<class InteractionVolume, class DataHandle>
590 void prepareDiffusionHandles_(InteractionVolume& iv,
591 DataHandle& handle,
592 bool forceUpdateAll)
593 {
594 for (unsigned int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
595 {
596 for (unsigned int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
597 {
598 // skip main component
599 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
600 if constexpr (!FluidSystem::isTracerFluidSystem())
601 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
602 continue;
603
604 // fill data in the handle
605 handle.diffusionHandle().setPhaseIndex(phaseIdx);
606 handle.diffusionHandle().setComponentIndex(compIdx);
607
608 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
609
610 // get instance of the interaction volume-local assembler
611 using Traits = typename InteractionVolume::Traits;
612 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
613 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
614
615 // maybe (re-)assemble matrices
616 if (forceUpdateAll || diffusionIsSolDependent)
617 {
618 // lambda to obtain diffusion coefficient
619 const auto getD = [&](const auto& volVars)
620 {
621 if constexpr (FluidSystem::isTracerFluidSystem())
622 return volVars.effectiveDiffusionCoefficient(0, 0, compIdx);
623 else
624 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
625 };
626
627 // Effective diffusion coefficients might be zero if saturation = 0.
628 // Compute epsilon to detect obsolete rows in the iv-local matrices during assembly
629 static const auto zeroD = getParamFromGroup<Scalar>(
630 problem().paramGroup(),
631 "Mpfa.ZeroEffectiveDiffusionCoefficientThreshold",
632 1e-16
633 );
634
635 // compute a representative transmissibility for this interaction volume, using
636 // the threshold for zero diffusion coefficients, and use this as epsilon
637 const auto& scv = fvGeometry().scv(iv.localScv(0).gridScvIndex());
638 const auto& scvf = fvGeometry().scvf(iv.localScvf(0).gridScvfIndex());
639 const auto& vv = elemVolVars()[scv];
640 const auto eps = Extrusion::area(scvf)*computeTpfaTransmissibility(
641 scvf, scv, zeroD, vv.extrusionFactor()
642 );
643
644 localAssembler.assembleMatrices(handle.diffusionHandle(), iv, getD, eps);
645 }
646
647 // assemble vector of mole fractions
648 auto getMassOrMoleFraction = [phaseIdx, compIdx] (const auto& volVars)
649 {
650 return (DiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged) ? volVars.massFraction(phaseIdx, compIdx) :
651 volVars.moleFraction(phaseIdx, compIdx);
652 };
653
654 localAssembler.assembleU(handle.diffusionHandle(), iv, getMassOrMoleFraction);
655 }
656 }
657 }
658
660 template<class InteractionVolume, class DataHandle>
661 void prepareHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
662 {
663 // get instance of the interaction volume-local assembler
664 using Traits = typename InteractionVolume::Traits;
665 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
666 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
667
668 // lambda to obtain the effective thermal conductivity
669 auto getLambda = [] (const auto& volVars) { return volVars.effectiveThermalConductivity(); };
670
671 // maybe (re-)assemble matrices
672 if (forceUpdateAll || heatConductionIsSolDependent)
673 localAssembler.assembleMatrices(handle.heatConductionHandle(), iv, getLambda);
674
675 // assemble vector of temperatures
676 auto getTemperature = [] (const auto& volVars) { return volVars.temperature(); };
677 localAssembler.assembleU(handle.heatConductionHandle(), iv, getTemperature);
678 }
679
680 const Problem* problemPtr_;
681 const FVElementGeometry* fvGeometryPtr_;
682 const ElementVolumeVariables* elemVolVarsPtr_;
683
684 // We store pointers to an inner and a boundary interaction volume.
685 // These are updated during the filling of the caches and the
686 // physics-related caches have access to them
687 PrimaryInteractionVolume* primaryIv_;
688 SecondaryInteractionVolume* secondaryIv_;
689
690 // pointer to the current interaction volume data handle
691 PrimaryDataHandle* primaryIvDataHandle_;
692 SecondaryDataHandle* secondaryIvDataHandle_;
693
694 // We do an interaction volume-wise filling of the caches
695 // While filling, we store a pointer to the current localScvf
696 // face data object of the IV so that the individual caches
697 // can access it and don't have to retrieve it again
698 const PrimaryLocalFaceData* primaryLocalFaceData_;
699 const SecondaryLocalFaceData* secondaryLocalFaceData_;
700};
701
702} // end namespace Dumux
703
704#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:47
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
Definition: propertysystem.hh:150
constexpr CCMpfa ccmpfa
Definition: method.hh:138
Definition: fluxvariablescachefiller.hh:39
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets the problem pointer.
Definition: fluxvariablescachefiller.hh:80
void fill(FluxVariablesCacheContainer &fluxVarsCacheContainer, FluxVariablesCache &scvfFluxVarsCache, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, bool forceUpdateAll=false)
function to fill the flux variables caches
Definition: fluxvariablescachefiller.hh:95
const PrimaryInteractionVolume & primaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:330
const SecondaryDataHandle & secondaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:342
const PrimaryDataHandle & primaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:338
const SecondaryLocalFaceData & secondaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:350
void fill(FluxVarsCacheStorage &fluxVarsCacheStorage, FluxVariablesCache &scvfFluxVarsCache, IVDataStorage &ivDataStorage, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, bool forceUpdateAll=false)
function to fill the flux variables caches
Definition: fluxvariablescachefiller.hh:247
const PrimaryLocalFaceData & primaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:346
const SecondaryInteractionVolume & secondaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:334
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets problem pointer.
Definition: fluxvariablescachefiller.hh:232
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...