3.2-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
30
34
35namespace Dumux {
36
37// forward declaration
38template<class TypeTag, DiscretizationMethod discMethod>
40
47template<class TypeTag>
49
51template<class TypeTag>
53{
57 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
58 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
59 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
60 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
61
62 using Element = typename GridView::template Codim<0>::Entity;
63
64 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
65 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
66 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
67
68 static constexpr bool advectionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
69 static constexpr bool diffusionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
70 static constexpr bool heatConductionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
71
72
73public:
74 static constexpr bool isSolDependent = (advectionEnabled && advectionIsSolDependent) ||
75 (diffusionEnabled && diffusionIsSolDependent) ||
76 (heatConductionEnabled && heatConductionIsSolDependent);
77
80 : problemPtr_(&problem) {}
81
93 template<class FluxVariablesCacheContainer, class FluxVariablesCache>
94 void fill(FluxVariablesCacheContainer& fluxVarsCacheContainer,
95 FluxVariablesCache& scvfFluxVarsCache,
96 const Element& element,
97 const FVElementGeometry& fvGeometry,
98 const ElementVolumeVariables& elemVolVars,
99 const SubControlVolumeFace& scvf,
100 bool forceUpdateAll = false)
101 {
102 // fill the physics-related quantities of the caches
103 if (forceUpdateAll)
104 {
105 if constexpr (advectionEnabled)
106 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
107 if constexpr (diffusionEnabled)
108 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
109 if constexpr (heatConductionEnabled)
110 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
111 }
112 else
113 {
114 if constexpr (advectionEnabled && advectionIsSolDependent)
115 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
116 if constexpr (diffusionEnabled && diffusionIsSolDependent)
117 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
118 if constexpr (heatConductionEnabled && heatConductionIsSolDependent)
119 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
120 }
121 }
122
123private:
124
125 const Problem& problem() const
126 { return *problemPtr_; }
127
129 template<class FluxVariablesCache>
130 void fillAdvection_(FluxVariablesCache& scvfFluxVarsCache,
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& elemVolVars,
134 const SubControlVolumeFace& scvf)
135 {
136 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
137 using AdvectionFiller = typename AdvectionType::Cache::Filler;
138
139 // forward to the filler for the advective quantities
140 AdvectionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
141 }
142
144 template<class FluxVariablesCache>
145 void fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
146 const Element& element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& elemVolVars,
149 const SubControlVolumeFace& scvf)
150 {
151 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
152 using DiffusionFiller = typename DiffusionType::Cache::Filler;
153 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
154
155 static constexpr int numPhases = ModelTraits::numFluidPhases();
156 static constexpr int numComponents = ModelTraits::numFluidComponents();
157
158 // forward to the filler of the diffusive quantities
159 if constexpr (FluidSystem::isTracerFluidSystem())
160 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
161 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
162 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *this);
163 else
164 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
165 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
166 if (compIdx != FluidSystem::getMainComponent(phaseIdx))
167 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *this);
168 }
169
171 template<class FluxVariablesCache>
172 void fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const SubControlVolumeFace& scvf)
177 {
178 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
179 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
180
181 // forward to the filler of the diffusive quantities
182 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
183 }
184
185 const Problem* problemPtr_;
186};
187
189template<class TypeTag>
191{
195 using Element = typename GridView::template Codim<0>::Entity;
196
198 using FVElementGeometry = typename GridGeometry::LocalView;
199 using MpfaHelper = typename GridGeometry::MpfaHelper;
200 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
201 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
202 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
203
204 using PrimaryInteractionVolume = GetPropType<TypeTag, Properties::PrimaryInteractionVolume>;
205 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle;
206 using PrimaryLocalFaceData = typename PrimaryInteractionVolume::Traits::LocalFaceData;
207 using SecondaryInteractionVolume = GetPropType<TypeTag, Properties::SecondaryInteractionVolume>;
208 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle;
209 using SecondaryLocalFaceData = typename SecondaryInteractionVolume::Traits::LocalFaceData;
210
211 static constexpr int dim = GridView::dimension;
212 static constexpr int dimWorld = GridView::dimensionworld;
213
214 static constexpr bool advectionEnabled = ModelTraits::enableAdvection();
215 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
216 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
217
218 static constexpr bool advectionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
219 static constexpr bool diffusionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
220 static constexpr bool heatConductionIsSolDependent = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
221
222public:
226 static constexpr bool isSolDependent = true;
227
230 : problemPtr_(&problem) {}
231
244 template<class FluxVarsCacheStorage, class FluxVariablesCache, class IVDataStorage>
245 void fill(FluxVarsCacheStorage& fluxVarsCacheStorage,
246 FluxVariablesCache& scvfFluxVarsCache,
247 IVDataStorage& ivDataStorage,
248 const Element& element,
249 const FVElementGeometry& fvGeometry,
250 const ElementVolumeVariables& elemVolVars,
251 const SubControlVolumeFace& scvf,
252 bool forceUpdateAll = false)
253 {
254 // Set pointers
255 elementPtr_ = &element;
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 Element& element() const { return *elementPtr_; }
357 const FVElementGeometry& fvGeometry() const { return *fvGeometryPtr_; }
358 const ElementVolumeVariables& elemVolVars() const { return *elemVolVarsPtr_; }
359
361 template<class FluxVariablesCache, class FluxVarsCacheStorage, class InteractionVolume>
362 void fillCachesInInteractionVolume_(FluxVarsCacheStorage& fluxVarsCacheStorage,
363 InteractionVolume& iv,
364 unsigned int ivIndexInContainer)
365 {
366 // determine if secondary interaction volumes are used here
367 static constexpr bool isSecondary = MpfaHelper::considerSecondaryIVs()
368 && std::is_same_v<InteractionVolume, SecondaryInteractionVolume>;
369
370 // First we upate data which are not dependent on the physical processes.
371 // We store pointers to the other flux var caches, so that we have to obtain
372 // this data only once and can use it again in the sub-cache fillers.
373 const auto numGlobalScvfs = iv.localFaceData().size();
374 std::vector<const SubControlVolumeFace*> ivScvfs(numGlobalScvfs);
375 std::vector<FluxVariablesCache*> ivFluxVarCaches(numGlobalScvfs);
376
377 unsigned int i = 0;
378 for (const auto& d : iv.localFaceData())
379 {
380 // obtain the scvf
381 const auto& scvfJ = fvGeometry().scvf(d.gridScvfIndex());
382 ivScvfs[i] = &scvfJ;
383 ivFluxVarCaches[i] = &fluxVarsCacheStorage[scvfJ];
384 ivFluxVarCaches[i]->setIvIndexInContainer(ivIndexInContainer);
385 ivFluxVarCaches[i]->setUpdateStatus(true);
386 ivFluxVarCaches[i]->setSecondaryIvUsage(isSecondary);
387 ivFluxVarCaches[i]->setIvLocalFaceIndex(d.ivLocalScvfIndex());
388 if (dim < dimWorld)
389 if (d.isOutsideFace())
390 ivFluxVarCaches[i]->setIndexInOutsideFaces(d.scvfLocalOutsideScvfIndex());
391 i++;
392 }
393
394 if constexpr (advectionEnabled)
395 fillAdvection_(iv, ivScvfs, ivFluxVarCaches);
396 if constexpr (diffusionEnabled)
397 fillDiffusion_(iv, ivScvfs, ivFluxVarCaches);
398 if constexpr (heatConductionEnabled)
399 fillHeatConduction_(iv, ivScvfs, ivFluxVarCaches);
400 }
401
403 template<class InteractionVolume, class FluxVariablesCache>
404 void fillAdvection_(InteractionVolume& iv,
405 const std::vector<const SubControlVolumeFace*>& ivScvfs,
406 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
407 {
408 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
409 using AdvectionFiller = typename AdvectionType::Cache::Filler;
410
411 // fill advection caches
412 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
413 {
414 // set pointer to current local face data object
415 // ifs are evaluated at compile time and are optimized away
416 if (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
417 {
418 // we cannot make a disctinction, thus we set both pointers
419 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
420 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
421 }
422 else if (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
423 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
424 else
425 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
426
427 // fill this scvfs cache
428 AdvectionFiller::fill(*ivFluxVarCaches[i],
429 problem(),
430 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
431 fvGeometry(),
432 elemVolVars(),
433 *ivScvfs[i],
434 *this);
435 }
436 }
437
439 template<class InteractionVolume, class FluxVariablesCache>
440 void fillDiffusion_(InteractionVolume& iv,
441 const std::vector<const SubControlVolumeFace*>& ivScvfs,
442 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
443 {
444 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
445 using DiffusionFiller = typename DiffusionType::Cache::Filler;
446
447 static constexpr int numPhases = ModelTraits::numFluidPhases();
448 static constexpr int numComponents = ModelTraits::numFluidComponents();
449
450 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
451 {
452 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
453 {
454 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
455 if constexpr (!FluidSystem::isTracerFluidSystem())
456 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
457 continue;
458
459 // fill diffusion caches
460 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
461 {
462 // set pointer to current local face data object
463 // ifs are evaluated at compile time and are optimized away
464 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
465 {
466 // we cannot make a disctinction, thus we set both pointers
467 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
468 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
469 }
470 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
471 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
472 else
473 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
474
475 // fill this scvfs cache
476 DiffusionFiller::fill(*ivFluxVarCaches[i],
477 phaseIdx,
478 compIdx,
479 problem(),
480 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
481 fvGeometry(),
482 elemVolVars(),
483 *ivScvfs[i],
484 *this);
485 }
486 }
487 }
488 }
489
491 template<class InteractionVolume, class FluxVariablesCache>
492 void fillHeatConduction_(InteractionVolume& iv,
493 const std::vector<const SubControlVolumeFace*>& ivScvfs,
494 const std::vector<FluxVariablesCache*>& ivFluxVarCaches)
495 {
496 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
497 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
498
499 // fill heat conduction caches
500 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
501 {
502 // set pointer to current local face data object
503 // ifs are evaluated at compile time and are optimized away
504 if constexpr (std::is_same_v<PrimaryInteractionVolume, SecondaryInteractionVolume>)
505 {
506 // we cannot make a disctinction, thus we set both pointers
507 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
508 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
509 }
510 else if constexpr (std::is_same_v<InteractionVolume, PrimaryInteractionVolume>)
511 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
512 else
513 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
514
515 // fill this scvfs cache
516 HeatConductionFiller::fill(*ivFluxVarCaches[i],
517 problem(),
518 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
519 fvGeometry(),
520 elemVolVars(),
521 *ivScvfs[i],
522 *this);
523 }
524 }
525
527 template< class InteractionVolume, class DataHandle>
528 void prepareDataHandle_([[maybe_unused]] InteractionVolume& iv, [[maybe_unused]] DataHandle& handle, [[maybe_unused]] bool forceUpdate)
529 {
530 // (maybe) solve system subject to intrinsic permeability
531 if constexpr (advectionEnabled)
532 {
533 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
534 if constexpr (AdvectionType::discMethod == DiscretizationMethod::ccmpfa)
535 prepareAdvectionHandle_(iv, handle, forceUpdate);
536 }
537
538 // (maybe) solve system subject to diffusion tensors
539 if constexpr (diffusionEnabled)
540 {
541 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
542 if constexpr (DiffusionType::discMethod == DiscretizationMethod::ccmpfa)
543 prepareDiffusionHandles_(iv, handle, forceUpdate);
544 }
545
546 // (maybe) solve system subject to thermal conductivity
547 if constexpr (heatConductionEnabled)
548 {
549 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
550 if constexpr (HeatConductionType::discMethod == DiscretizationMethod::ccmpfa)
551 prepareHeatConductionHandle_(iv, handle, forceUpdate);
552 }
553 }
554
556 template<class InteractionVolume, class DataHandle>
557 void prepareAdvectionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
558 {
559 // get instance of the interaction volume-local assembler
560 using Traits = typename InteractionVolume::Traits;
561 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
562 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
563
564 // lambda to obtain the permeability tensor
565 auto getK = [] (const auto& volVars) { return volVars.permeability(); };
566
567 // Assemble T only if permeability is sol-dependent or if update is forced
568 if (forceUpdateAll || advectionIsSolDependent)
569 localAssembler.assembleMatrices(handle.advectionHandle(), iv, getK);
570
571 // assemble pressure vectors
572 for (unsigned int pIdx = 0; pIdx < ModelTraits::numFluidPhases(); ++pIdx)
573 {
574 // set context in handle
575 handle.advectionHandle().setPhaseIndex(pIdx);
576
577 // maybe (re-)assemble gravity contribution vector
578 auto getRho = [pIdx] (const auto& volVars) { return volVars.density(pIdx); };
579 static const bool enableGravity = getParamFromGroup<bool>(problem().paramGroup(), "Problem.EnableGravity");
580 if (enableGravity)
581 localAssembler.assembleGravity(handle.advectionHandle(), iv, getRho);
582
583 // reassemble pressure vector
584 auto getPressure = [pIdx] (const auto& volVars) { return volVars.pressure(pIdx); };
585 localAssembler.assembleU(handle.advectionHandle(), iv, getPressure);
586 }
587 }
588
590 template<class InteractionVolume, class DataHandle>
591 void prepareDiffusionHandles_(InteractionVolume& iv,
592 DataHandle& handle,
593 bool forceUpdateAll)
594 {
595 for (unsigned int phaseIdx = 0; phaseIdx < ModelTraits::numFluidPhases(); ++phaseIdx)
596 {
597 for (unsigned int compIdx = 0; compIdx < ModelTraits::numFluidComponents(); ++compIdx)
598 {
599 // skip main component
600 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
601 if constexpr (!FluidSystem::isTracerFluidSystem())
602 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
603 continue;
604
605 // fill data in the handle
606 handle.diffusionHandle().setPhaseIndex(phaseIdx);
607 handle.diffusionHandle().setComponentIndex(compIdx);
608
609 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
610 using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
611
612 // get instance of the interaction volume-local assembler
613 using Traits = typename InteractionVolume::Traits;
614 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
615 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
616
617 // maybe (re-)assemble matrices
618 if (forceUpdateAll || diffusionIsSolDependent)
619 {
620 // lambda to obtain diffusion coefficient
621 const auto getD = [&](const auto& volVars)
622 {
623 if constexpr (FluidSystem::isTracerFluidSystem())
624 return Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(volVars, 0, 0, compIdx);
625 else
626 return Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(volVars, phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
627 };
628
629 // Effective diffusion coefficients might get zero if saturation = 0.
630 // Compute epsilon to detect obsolete rows in the iv-local matrices during assembly
631 const auto& scv = *scvs(fvGeometry()).begin();
632 const auto& scvf = *scvfs(fvGeometry()).begin();
633 const auto& vv = elemVolVars()[scv];
634 const auto D = [&] ()
635 {
636 // diffusion coefficients below 1e-20 are treated as zeroes!!
637 using std::max;
638 if constexpr (!FluidSystem::isTracerFluidSystem())
639 return max(1e-20, vv.diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx));
640 else
641 return max(1e-20, vv.diffusionCoefficient(0, 0, compIdx));
642 } ();
643
644 auto eps = 1e-7*computeTpfaTransmissibility(scvf, scv, D, vv.extrusionFactor())*scvf.area();
645 localAssembler.assembleMatrices(handle.diffusionHandle(), iv, getD, eps);
646 }
647
648 // assemble vector of mole fractions
649 auto getMassOrMoleFraction = [phaseIdx, compIdx] (const auto& volVars)
650 {
651 return (DiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged) ? volVars.massFraction(phaseIdx, compIdx) :
652 volVars.moleFraction(phaseIdx, compIdx);
653 };
654
655 localAssembler.assembleU(handle.diffusionHandle(), iv, getMassOrMoleFraction);
656 }
657 }
658 }
659
661 template<class InteractionVolume, class DataHandle>
662 void prepareHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
663 {
664 // get instance of the interaction volume-local assembler
665 using Traits = typename InteractionVolume::Traits;
666 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
667 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
668
669 // lambda to obtain the effective thermal conductivity
670 auto getLambda = [] (const auto& volVars) { return volVars.effectiveThermalConductivity(); };
671
672 // maybe (re-)assemble matrices
673 if (forceUpdateAll || heatConductionIsSolDependent)
674 localAssembler.assembleMatrices(handle.heatConductionHandle(), iv, getLambda);
675
676 // assemble vector of temperatures
677 auto getTemperature = [] (const auto& volVars) { return volVars.temperature(); };
678 localAssembler.assembleU(handle.heatConductionHandle(), iv, getTemperature);
679 }
680
681 const Problem* problemPtr_;
682 const Element* elementPtr_;
683 const FVElementGeometry* fvGeometryPtr_;
684 const ElementVolumeVariables* elemVolVarsPtr_;
685
686 // We store pointers to an inner and a boundary interaction volume.
687 // These are updated during the filling of the caches and the
688 // physics-related caches have access to them
689 PrimaryInteractionVolume* primaryIv_;
690 SecondaryInteractionVolume* secondaryIv_;
691
692 // pointer to the current interaction volume data handle
693 PrimaryDataHandle* primaryIvDataHandle_;
694 SecondaryDataHandle* secondaryIvDataHandle_;
695
696 // We do an interaction volume-wise filling of the caches
697 // While filling, we store a pointer to the current localScvf
698 // face data object of the IV so that the individual caches
699 // can access it and don't have to retrieve it again
700 const PrimaryLocalFaceData* primaryLocalFaceData_;
701 const SecondaryLocalFaceData* secondaryLocalFaceData_;
702};
703
704} // end namespace Dumux
705
706#endif
Helpers for deprecation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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
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 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:79
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:94
const PrimaryLocalFaceData & primaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:346
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:245
const PrimaryInteractionVolume & primaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:330
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets problem pointer.
Definition: fluxvariablescachefiller.hh:229
const PrimaryDataHandle & primaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:338
const SecondaryDataHandle & secondaryIvDataHandle() const
returns the stored data handle pointer
Definition: fluxvariablescachefiller.hh:342
const SecondaryLocalFaceData & secondaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: fluxvariablescachefiller.hh:350
const SecondaryInteractionVolume & secondaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: fluxvariablescachefiller.hh:334
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...