3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/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
33
34namespace Dumux {
35
36// forward declaration
37template<class TypeTag, DiscretizationMethod discMethod>
39
46template<class TypeTag>
48
50template<class TypeTag>
52{
56 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
57 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
59 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
60
61 using Element = typename GridView::template Codim<0>::Entity;
62
63 static constexpr bool doAdvection = ModelTraits::enableAdvection();
64 static constexpr bool doDiffusion = ModelTraits::enableMolecularDiffusion();
65 static constexpr bool doHeatConduction = ModelTraits::enableEnergyBalance();
66
67 static constexpr bool soldependentAdvection = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
68 static constexpr bool soldependentDiffusion = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
69 static constexpr bool soldependentHeatConduction = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
70
71
72public:
73 static constexpr bool isSolDependent = (doAdvection && soldependentAdvection) ||
74 (doDiffusion && soldependentDiffusion) ||
75 (doHeatConduction && soldependentHeatConduction);
76
79 : problemPtr_(&problem) {}
80
92 template<class FluxVariablesCacheContainer, class FluxVariablesCache>
93 void fill(FluxVariablesCacheContainer& fluxVarsCacheContainer,
94 FluxVariablesCache& scvfFluxVarsCache,
95 const Element& element,
96 const FVElementGeometry& fvGeometry,
97 const ElementVolumeVariables& elemVolVars,
98 const SubControlVolumeFace& scvf,
99 bool forceUpdateAll = false)
100 {
101 // fill the physics-related quantities of the caches
102 if (forceUpdateAll)
103 {
104 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
105 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
106 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
107 }
108 else
109 {
110 if (doAdvection && soldependentAdvection)
111 fillAdvection_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
112 if (doDiffusion && soldependentDiffusion)
113 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
114 if (doHeatConduction && soldependentHeatConduction)
115 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
116 }
117 }
118
119private:
120
121 const Problem& problem() const
122 { return *problemPtr_; }
123
125 template<class FluxVariablesCache, bool advectionEnabled = doAdvection>
126 typename std::enable_if<advectionEnabled>::type
127 fillAdvection_(FluxVariablesCache& scvfFluxVarsCache,
128 const Element& element,
129 const FVElementGeometry& fvGeometry,
130 const ElementVolumeVariables& elemVolVars,
131 const SubControlVolumeFace& scvf)
132 {
133 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
134 using AdvectionFiller = typename AdvectionType::Cache::Filler;
135
136 // forward to the filler for the advective quantities
137 AdvectionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
138 }
139
141 template<class FluxVariablesCache, bool advectionEnabled = doAdvection>
142 typename std::enable_if<!advectionEnabled>::type
143 fillAdvection_(FluxVariablesCache& scvfFluxVarsCache,
144 const Element& element,
145 const FVElementGeometry& fvGeometry,
146 const ElementVolumeVariables& elemVolVars,
147 const SubControlVolumeFace& scvf)
148 {}
149
151 template<class FluxVariablesCache, bool diffusionEnabled = doDiffusion>
152 typename std::enable_if<diffusionEnabled>::type
153 fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
154 const Element& element,
155 const FVElementGeometry& fvGeometry,
156 const ElementVolumeVariables& elemVolVars,
157 const SubControlVolumeFace& scvf)
158 {
159 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
160 using DiffusionFiller = typename DiffusionType::Cache::Filler;
161 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
162
163 static constexpr int numPhases = ModelTraits::numFluidPhases();
164 static constexpr int numComponents = ModelTraits::numFluidComponents();
165
166 // forward to the filler of the diffusive quantities
167 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
168 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
169 if (compIdx != FluidSystem::getMainComponent(phaseIdx))
170 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *this);
171 }
172
174 template<class FluxVariablesCache, bool diffusionEnabled = doDiffusion>
175 typename std::enable_if<!diffusionEnabled>::type
176 fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
177 const Element& element,
178 const FVElementGeometry& fvGeometry,
179 const ElementVolumeVariables& elemVolVars,
180 const SubControlVolumeFace& scvf)
181 {}
182
184 template<class FluxVariablesCache, bool heatConductionEnabled = doHeatConduction>
185 typename std::enable_if<heatConductionEnabled>::type
186 fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
187 const Element& element,
188 const FVElementGeometry& fvGeometry,
189 const ElementVolumeVariables& elemVolVars,
190 const SubControlVolumeFace& scvf)
191 {
192 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
193 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
194
195 // forward to the filler of the diffusive quantities
196 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
197 }
198
200 template<class FluxVariablesCache, bool heatConductionEnabled = doHeatConduction>
201 typename std::enable_if<!heatConductionEnabled>::type
202 fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
203 const Element& element,
204 const FVElementGeometry& fvGeometry,
205 const ElementVolumeVariables& elemVolVars,
206 const SubControlVolumeFace& scvf)
207 {}
208
209 const Problem* problemPtr_;
210};
211
213template<class TypeTag>
215{
219 using Element = typename GridView::template Codim<0>::Entity;
220
222 using FVElementGeometry = typename GridGeometry::LocalView;
223 using MpfaHelper = typename GridGeometry::MpfaHelper;
224 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
225 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
226 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
227
228 using PrimaryInteractionVolume = GetPropType<TypeTag, Properties::PrimaryInteractionVolume>;
229 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle;
230 using PrimaryLocalFaceData = typename PrimaryInteractionVolume::Traits::LocalFaceData;
231 using SecondaryInteractionVolume = GetPropType<TypeTag, Properties::SecondaryInteractionVolume>;
232 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle;
233 using SecondaryLocalFaceData = typename SecondaryInteractionVolume::Traits::LocalFaceData;
234
235 static constexpr int dim = GridView::dimension;
236 static constexpr int dimWorld = GridView::dimensionworld;
237
238 static constexpr bool doAdvection = ModelTraits::enableAdvection();
239 static constexpr bool doDiffusion = ModelTraits::enableMolecularDiffusion();
240 static constexpr bool doHeatConduction = ModelTraits::enableEnergyBalance();
241
242 static constexpr bool soldependentAdvection = getPropValue<TypeTag, Properties::SolutionDependentAdvection>();
243 static constexpr bool soldependentDiffusion = getPropValue<TypeTag, Properties::SolutionDependentMolecularDiffusion>();
244 static constexpr bool soldependentHeatConduction = getPropValue<TypeTag, Properties::SolutionDependentHeatConduction>();
245
246public:
250 static constexpr bool isSolDependent = true;
251
254 : problemPtr_(&problem) {}
255
268 template<class FluxVarsCacheStorage, class FluxVariablesCache, class IVDataStorage>
269 void fill(FluxVarsCacheStorage& fluxVarsCacheStorage,
270 FluxVariablesCache& scvfFluxVarsCache,
271 IVDataStorage& ivDataStorage,
272 const Element& element,
273 const FVElementGeometry& fvGeometry,
274 const ElementVolumeVariables& elemVolVars,
275 const SubControlVolumeFace& scvf,
276 bool forceUpdateAll = false)
277 {
278 // Set pointers
279 elementPtr_ = &element;
280 fvGeometryPtr_ = &fvGeometry;
281 elemVolVarsPtr_ = &elemVolVars;
282
283 // prepare interaction volume and fill caches of all the scvfs connected to it
284 const auto& gridGeometry = fvGeometry.gridGeometry();
285 if (gridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
286 {
287 if (forceUpdateAll)
288 {
289 // the local index of the interaction volume to be created in its container
290 const auto ivIndexInContainer = ivDataStorage.secondaryInteractionVolumes.size();
291
292 // prepare the locally cached boundary interaction volume
293 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().secondaryIndexSet(scvf);
294 ivDataStorage.secondaryInteractionVolumes.emplace_back();
295 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes.back();
296 secondaryIv_->bind(indexSet, problem(), fvGeometry);
297
298 // prepare the corresponding data handle
299 ivDataStorage.secondaryDataHandles.emplace_back();
300 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles.back();
301
302 // fill the caches for all the scvfs in the interaction volume
303 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, *secondaryIvDataHandle_, ivIndexInContainer, true);
304 }
305 else
306 {
307 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
308 secondaryIv_ = &ivDataStorage.secondaryInteractionVolumes[ivIndexInContainer];
309 secondaryIvDataHandle_ = &ivDataStorage.secondaryDataHandles[ivIndexInContainer];
310
311 // fill the caches for all the scvfs in the interaction volume
312 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *secondaryIv_, *secondaryIvDataHandle_, ivIndexInContainer);
313 }
314 }
315 else
316 {
317 if (forceUpdateAll)
318 {
319 // the local index of the interaction volume to be created in its container
320 const auto ivIndexInContainer = ivDataStorage.primaryInteractionVolumes.size();
321
322 // prepare the locally cached boundary interaction volume
323 const auto& indexSet = gridGeometry.gridInteractionVolumeIndexSets().primaryIndexSet(scvf);
324 ivDataStorage.primaryInteractionVolumes.emplace_back();
325 primaryIv_ = &ivDataStorage.primaryInteractionVolumes.back();
326 primaryIv_->bind(indexSet, problem(), fvGeometry);
327
328 // prepare the corresponding data handle
329 ivDataStorage.primaryDataHandles.emplace_back();
330 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles.back();
331
332 // fill the caches for all the scvfs in the interaction volume
333 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, *primaryIvDataHandle_, ivIndexInContainer, true);
334 }
335 else
336 {
337 const auto ivIndexInContainer = scvfFluxVarsCache.ivIndexInContainer();
338 primaryIv_ = &ivDataStorage.primaryInteractionVolumes[ivIndexInContainer];
339 primaryIvDataHandle_ = &ivDataStorage.primaryDataHandles[ivIndexInContainer];
340
341 // fill the caches for all the scvfs in the interaction volume
342 fillCachesInInteractionVolume_<FluxVariablesCache>(fluxVarsCacheStorage, *primaryIv_, *primaryIvDataHandle_, ivIndexInContainer);
343 }
344 }
345 }
346
348 const PrimaryInteractionVolume& primaryInteractionVolume() const
349 { return *primaryIv_; }
350
352 const SecondaryInteractionVolume& secondaryInteractionVolume() const
353 { return *secondaryIv_; }
354
356 const PrimaryDataHandle& primaryIvDataHandle() const
357 { return *primaryIvDataHandle_; }
358
360 const SecondaryDataHandle& secondaryIvDataHandle() const
361 { return *secondaryIvDataHandle_; }
362
364 const PrimaryLocalFaceData& primaryIvLocalFaceData() const
365 { return *primaryLocalFaceData_; }
366
368 const SecondaryLocalFaceData& secondaryIvLocalFaceData() const
369 { return *secondaryLocalFaceData_; }
370
371private:
372
373 const Problem& problem() const { return *problemPtr_; }
374 const Element& element() const { return *elementPtr_; }
375 const FVElementGeometry& fvGeometry() const { return *fvGeometryPtr_; }
376 const ElementVolumeVariables& elemVolVars() const { return *elemVolVarsPtr_; }
377
379 template<class FluxVariablesCache, class FluxVarsCacheStorage, class InteractionVolume, class DataHandle>
380 void fillCachesInInteractionVolume_(FluxVarsCacheStorage& fluxVarsCacheStorage,
381 InteractionVolume& iv,
382 DataHandle& handle,
383 unsigned int ivIndexInContainer,
384 bool forceUpdateAll = false)
385 {
386 // determine if secondary interaction volumes are used here
387 static constexpr bool isSecondary = MpfaHelper::considerSecondaryIVs()
388 && std::is_same<InteractionVolume, SecondaryInteractionVolume>::value;
389
390 // First we upate data which are not dependent on the physical processes.
391 // We store pointers to the other flux var caches, so that we have to obtain
392 // this data only once and can use it again in the sub-cache fillers.
393 const auto numGlobalScvfs = iv.localFaceData().size();
394 std::vector<const SubControlVolumeFace*> ivScvfs(numGlobalScvfs);
395 std::vector<FluxVariablesCache*> ivFluxVarCaches(numGlobalScvfs);
396
397 unsigned int i = 0;
398 for (const auto& d : iv.localFaceData())
399 {
400 // obtain the scvf
401 const auto& scvfJ = fvGeometry().scvf(d.gridScvfIndex());
402 ivScvfs[i] = &scvfJ;
403 ivFluxVarCaches[i] = &fluxVarsCacheStorage[scvfJ];
404 ivFluxVarCaches[i]->setIvIndexInContainer(ivIndexInContainer);
405 ivFluxVarCaches[i]->setUpdateStatus(true);
406 ivFluxVarCaches[i]->setSecondaryIvUsage(isSecondary);
407 ivFluxVarCaches[i]->setIvLocalFaceIndex(d.ivLocalScvfIndex());
408 if (dim < dimWorld)
409 if (d.isOutsideFace())
410 ivFluxVarCaches[i]->setIndexInOutsideFaces(d.scvfLocalOutsideScvfIndex());
411 i++;
412 }
413
414 fillAdvection_(iv, handle, ivScvfs, ivFluxVarCaches, forceUpdateAll);
415 fillDiffusion_(iv, handle, ivScvfs, ivFluxVarCaches, forceUpdateAll);
416 fillHeatConduction_(iv, handle, ivScvfs, ivFluxVarCaches, forceUpdateAll);
417 }
418
420 template< class InteractionVolume,
421 class DataHandle,
422 class FluxVariablesCache,
423 bool enableAdvection = doAdvection,
424 typename std::enable_if_t<enableAdvection, int> = 0 >
425 void fillAdvection_(InteractionVolume& iv,
426 DataHandle& handle,
427 const std::vector<const SubControlVolumeFace*>& ivScvfs,
428 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
429 bool forceUpdateAll = false)
430 {
431 using AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>;
432 using AdvectionFiller = typename AdvectionType::Cache::Filler;
433
434 // fill data in the handle
435 fillAdvectionHandle_(iv, handle, forceUpdateAll);
436
437 // fill advection caches
438 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
439 {
440 // set pointer to current local face data object
441 // ifs are evaluated at compile time and are optimized away
442 if (std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value)
443 {
444 // we cannot make a disctinction, thus we set both pointers
445 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
446 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
447 }
448 else if (std::is_same<InteractionVolume, PrimaryInteractionVolume>::value)
449 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
450 else
451 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
452
453 // fill this scvfs cache
454 AdvectionFiller::fill(*ivFluxVarCaches[i],
455 problem(),
456 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
457 fvGeometry(),
458 elemVolVars(),
459 *ivScvfs[i],
460 *this);
461 }
462 }
463
465 template< class InteractionVolume,
466 class DataHandle,
467 class FluxVariablesCache,
468 bool enableAdvection = doAdvection,
469 typename std::enable_if_t<!enableAdvection, int> = 0 >
470 void fillAdvection_(InteractionVolume& iv,
471 DataHandle& handle,
472 const std::vector<const SubControlVolumeFace*>& ivScvfs,
473 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
474 bool forceUpdateAll = false)
475 {}
476
478 template< class InteractionVolume,
479 class DataHandle,
480 class FluxVariablesCache,
481 bool enableDiffusion = doDiffusion,
482 typename std::enable_if_t<enableDiffusion, int> = 0 >
483 void fillDiffusion_(InteractionVolume& iv,
484 DataHandle& handle,
485 const std::vector<const SubControlVolumeFace*>& ivScvfs,
486 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
487 bool forceUpdateAll = false)
488 {
489 using DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>;
490 using DiffusionFiller = typename DiffusionType::Cache::Filler;
491
492 static constexpr int numPhases = ModelTraits::numFluidPhases();
493 static constexpr int numComponents = ModelTraits::numFluidComponents();
494
495 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
496 {
497 handle.diffusionHandle().setPhaseIndex(phaseIdx);
498 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
499 {
500 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
501 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
502 continue;
503
504 // fill data in the handle
505 handle.diffusionHandle().setComponentIndex(compIdx);
506 fillDiffusionHandle_(iv, handle, forceUpdateAll, phaseIdx, compIdx);
507
508 // fill diffusion caches
509 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
510 {
511 // set pointer to current local face data object
512 // ifs are evaluated at compile time and are optimized away
513 if (std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value)
514 {
515 // we cannot make a disctinction, thus we set both pointers
516 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
517 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
518 }
519 else if (std::is_same<InteractionVolume, PrimaryInteractionVolume>::value)
520 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
521 else
522 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
523
524 // fill this scvfs cache
525 DiffusionFiller::fill(*ivFluxVarCaches[i],
526 phaseIdx,
527 compIdx,
528 problem(),
529 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
530 fvGeometry(),
531 elemVolVars(),
532 *ivScvfs[i],
533 *this);
534 }
535 }
536 }
537 }
538
540 template< class InteractionVolume,
541 class DataHandle,
542 class FluxVariablesCache,
543 bool enableDiffusion = doDiffusion,
544 typename std::enable_if_t<!enableDiffusion, int> = 0 >
545 void fillDiffusion_(InteractionVolume& iv,
546 DataHandle& handle,
547 const std::vector<const SubControlVolumeFace*>& ivScvfs,
548 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
549 bool forceUpdateAll = false)
550 {}
551
553 template< class InteractionVolume,
554 class DataHandle,
555 class FluxVariablesCache,
556 bool enableHeatConduction = doHeatConduction,
557 typename std::enable_if_t<enableHeatConduction, int> = 0 >
558 void fillHeatConduction_(InteractionVolume& iv,
559 DataHandle& handle,
560 const std::vector<const SubControlVolumeFace*>& ivScvfs,
561 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
562 bool forceUpdateAll = false)
563 {
564 using HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>;
565 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
566
567 // prepare data in handle
568 fillHeatConductionHandle_(iv, handle, forceUpdateAll);
569
570 // fill heat conduction caches
571 for (unsigned int i = 0; i < iv.localFaceData().size(); ++i)
572 {
573 // set pointer to current local face data object
574 // ifs are evaluated at compile time and are optimized away
575 if (std::is_same<PrimaryInteractionVolume, SecondaryInteractionVolume>::value)
576 {
577 // we cannot make a disctinction, thus we set both pointers
578 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
579 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
580 }
581 else if (std::is_same<InteractionVolume, PrimaryInteractionVolume>::value)
582 primaryLocalFaceData_ = &(iv.localFaceData()[i]);
583 else
584 secondaryLocalFaceData_ = &(iv.localFaceData()[i]);
585
586 // fill this scvfs cache
587 HeatConductionFiller::fill(*ivFluxVarCaches[i],
588 problem(),
589 iv.element(iv.localFaceData()[i].ivLocalInsideScvIndex()),
590 fvGeometry(),
591 elemVolVars(),
592 *ivScvfs[i],
593 *this);
594 }
595 }
596
598 template< class InteractionVolume,
599 class DataHandle,
600 class FluxVariablesCache,
601 bool enableHeatConduction = doHeatConduction,
602 typename std::enable_if_t<!enableHeatConduction, int> = 0 >
603 void fillHeatConduction_(InteractionVolume& iv,
604 DataHandle& handle,
605 const std::vector<const SubControlVolumeFace*>& ivScvfs,
606 const std::vector<FluxVariablesCache*>& ivFluxVarCaches,
607 bool forceUpdateAll = false)
608 {}
609
611 template< class InteractionVolume,
612 class DataHandle,
613 class AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>,
614 typename std::enable_if_t<AdvectionType::discMethod == DiscretizationMethod::ccmpfa, int> = 0 >
615 void fillAdvectionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
616 {
617 using LambdaFactory = TensorLambdaFactory<DiscretizationMethod::ccmpfa>;
618
619 // get instance of the interaction volume-local assembler
620 using Traits = typename InteractionVolume::Traits;
621 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
622 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
623
624 // Assemble T only if permeability is sol-dependent or if update is forced
625 if (forceUpdateAll || soldependentAdvection)
626 localAssembler.assembleMatrices(handle.advectionHandle(), iv, LambdaFactory::getAdvectionLambda());
627
628 // assemble pressure vectors
629 for (unsigned int pIdx = 0; pIdx < ModelTraits::numFluidPhases(); ++pIdx)
630 {
631 // set context in handle
632 handle.advectionHandle().setPhaseIndex(pIdx);
633
634 // maybe (re-)assemble gravity contribution vector
635 auto getRho = [pIdx] (const auto& volVars) { return volVars.density(pIdx); };
636 static const bool enableGravity = getParamFromGroup<bool>(problem().paramGroup(), "Problem.EnableGravity");
637 if (enableGravity)
638 localAssembler.assembleGravity(handle.advectionHandle(), iv, getRho);
639
640 // reassemble pressure vector
641 auto getPressure = [pIdx] (const auto& volVars) { return volVars.pressure(pIdx); };
642 localAssembler.assembleU(handle.advectionHandle(), iv, getPressure);
643 }
644 }
645
647 template< class InteractionVolume,
648 class DataHandle,
649 class DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>,
650 typename std::enable_if_t<DiffusionType::discMethod == DiscretizationMethod::ccmpfa, int> = 0 >
651 void fillDiffusionHandle_(InteractionVolume& iv,
652 DataHandle& handle,
653 bool forceUpdateAll,
654 int phaseIdx, int compIdx)
655 {
656 using LambdaFactory = TensorLambdaFactory<DiscretizationMethod::ccmpfa>;
657
658 // get instance of the interaction volume-local assembler
659 using Traits = typename InteractionVolume::Traits;
660 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
661 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
662
663 // maybe (re-)assemble matrices
664 if (forceUpdateAll || soldependentDiffusion)
665 localAssembler.assembleMatrices(handle.diffusionHandle(),
666 iv,
667 LambdaFactory::getDiffusionLambda(phaseIdx, compIdx));
668
669 // assemble vector of mole fractions
670 auto getMassOrMoleFraction = [phaseIdx, compIdx] (const auto& volVars) { return (DiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged) ?volVars.massFraction(phaseIdx, compIdx) : volVars.moleFraction(phaseIdx, compIdx); };
671 localAssembler.assembleU(handle.diffusionHandle(), iv, getMassOrMoleFraction);
672 }
673
675 template< class InteractionVolume,
676 class DataHandle,
677 class HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>,
678 typename std::enable_if_t<HeatConductionType::discMethod == DiscretizationMethod::ccmpfa, int> = 0 >
679 void fillHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll)
680 {
681 using LambdaFactory = TensorLambdaFactory<DiscretizationMethod::ccmpfa>;
682 using ThermCondModel = GetPropType<TypeTag, Properties::ThermalConductivityModel>;
683
684 // get instance of the interaction volume-local assembler
685 using Traits = typename InteractionVolume::Traits;
686 using IvLocalAssembler = typename Traits::template LocalAssembler<Problem, FVElementGeometry, ElementVolumeVariables>;
687 IvLocalAssembler localAssembler(problem(), fvGeometry(), elemVolVars());
688
689 // maybe (re-)assemble matrices
690 if (forceUpdateAll || soldependentHeatConduction)
691 localAssembler.assembleMatrices(handle.heatConductionHandle(),
692 iv,
693 LambdaFactory::template getHeatConductionLambda<ThermCondModel>());
694
695 // assemble vector of temperatures
696 auto getTemperature = [] (const auto& volVars) { return volVars.temperature(); };
697 localAssembler.assembleU(handle.heatConductionHandle(), iv, getTemperature);
698 }
699
701 template< class InteractionVolume,
702 class DataHandle,
703 class AdvectionType = GetPropType<TypeTag, Properties::AdvectionType>,
704 typename std::enable_if_t<AdvectionType::discMethod != DiscretizationMethod::ccmpfa, int> = 0 >
705 void fillAdvectionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll) {}
706
708 template< class InteractionVolume,
709 class DataHandle,
710 class DiffusionType = GetPropType<TypeTag, Properties::MolecularDiffusionType>,
711 typename std::enable_if_t<DiffusionType::discMethod != DiscretizationMethod::ccmpfa, int> = 0 >
712 void fillDiffusionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll, int phaseIdx, int compIdx) {}
713
715 template< class InteractionVolume,
716 class DataHandle,
717 class HeatConductionType = GetPropType<TypeTag, Properties::HeatConductionType>,
718 typename std::enable_if_t<HeatConductionType::discMethod != DiscretizationMethod::ccmpfa, int> = 0 >
719 void fillHeatConductionHandle_(InteractionVolume& iv, DataHandle& handle, bool forceUpdateAll) {}
720
721 const Problem* problemPtr_;
722 const Element* elementPtr_;
723 const FVElementGeometry* fvGeometryPtr_;
724 const ElementVolumeVariables* elemVolVarsPtr_;
725
726 // We store pointers to an inner and a boundary interaction volume.
727 // These are updated during the filling of the caches and the
728 // physics-related caches have access to them
729 PrimaryInteractionVolume* primaryIv_;
730 SecondaryInteractionVolume* secondaryIv_;
731
732 // pointer to the current interaction volume data handle
733 PrimaryDataHandle* primaryIvDataHandle_;
734 SecondaryDataHandle* secondaryIvDataHandle_;
735
736 // We do an interaction volume-wise filling of the caches
737 // While filling, we store a pointer to the current localScvf
738 // face data object of the IV so that the individual caches
739 // can access it and don't have to retrieve it again
740 const PrimaryLocalFaceData* primaryLocalFaceData_;
741 const SecondaryLocalFaceData* secondaryLocalFaceData_;
742};
743
744} // end namespace Dumux
745
746#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper class to be used to obtain lambda functions for the tensors involved in the laws that describe...
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Definition: porousmediumflow/fluxvariablescachefiller.hh:38
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets the problem pointer.
Definition: porousmediumflow/fluxvariablescachefiller.hh:78
void fill(FluxVariablesCacheContainer &fluxVarsCacheContainer, FluxVariablesCache &scvfFluxVarsCache, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, bool forceUpdateAll=false)
function to fill the flux variables caches
Definition: porousmediumflow/fluxvariablescachefiller.hh:93
const PrimaryLocalFaceData & primaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: porousmediumflow/fluxvariablescachefiller.hh:364
void fill(FluxVarsCacheStorage &fluxVarsCacheStorage, FluxVariablesCache &scvfFluxVarsCache, IVDataStorage &ivDataStorage, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, bool forceUpdateAll=false)
function to fill the flux variables caches
Definition: porousmediumflow/fluxvariablescachefiller.hh:269
const PrimaryInteractionVolume & primaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: porousmediumflow/fluxvariablescachefiller.hh:348
PorousMediumFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets problem pointer.
Definition: porousmediumflow/fluxvariablescachefiller.hh:253
const PrimaryDataHandle & primaryIvDataHandle() const
returns the stored data handle pointer
Definition: porousmediumflow/fluxvariablescachefiller.hh:356
const SecondaryDataHandle & secondaryIvDataHandle() const
returns the stored data handle pointer
Definition: porousmediumflow/fluxvariablescachefiller.hh:360
const SecondaryLocalFaceData & secondaryIvLocalFaceData() const
returns the currently stored iv-local face data object
Definition: porousmediumflow/fluxvariablescachefiller.hh:368
const SecondaryInteractionVolume & secondaryInteractionVolume() const
returns the stored interaction volume pointer
Definition: porousmediumflow/fluxvariablescachefiller.hh:352
Declares all properties used in Dumux.