version 3.10-dev
scalarfluxvariablescachefiller.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_FREEFLOW_SCALAR_FLUXVARIABLESCACHE_FILLER_HH
13#define DUMUX_FREEFLOW_SCALAR_FLUXVARIABLESCACHE_FILLER_HH
14
17
22
23namespace Dumux {
24
25// forward declaration
26template<class Problem, class ModelTraits, bool diffusionIsSolDependent, bool heatConductionIsSolDependent, class DiscretizationMethod>
28
35template<class Problem, class ModelTraits, bool diffusionIsSolDependent, bool heatConductionIsSolDependent>
37
39template<class Problem, class ModelTraits, bool diffusionIsSolDependent, bool heatConductionIsSolDependent>
40class FreeFlowScalarFluxVariablesCacheFillerImplementation<Problem, ModelTraits, diffusionIsSolDependent, heatConductionIsSolDependent, DiscretizationMethods::CCTpfa>
41{
42 using GridGeometry = typename ProblemTraits<Problem>::GridGeometry;
43 using GridView = typename GridGeometry::GridView;
44 using FVElementGeometry = typename GridGeometry::LocalView;
45 using SubControlVolume = typename GridGeometry::SubControlVolume;
46 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
47
48 using Element = typename GridView::template Codim<0>::Entity;
49
50 static constexpr bool diffusionEnabled = ModelTraits::enableMolecularDiffusion();
51 static constexpr bool heatConductionEnabled = ModelTraits::enableEnergyBalance();
52
53public:
54 static constexpr bool isSolDependent = (diffusionEnabled && diffusionIsSolDependent) ||
55 (heatConductionEnabled && heatConductionIsSolDependent);
56
59 : problemPtr_(&problem) {}
60
72 template<class FluxVariablesCacheContainer, class FluxVariablesCache, class ElementVolumeVariables>
73 void fill(FluxVariablesCacheContainer& fluxVarsCacheContainer,
74 FluxVariablesCache& scvfFluxVarsCache,
75 const Element& element,
76 const FVElementGeometry& fvGeometry,
77 const ElementVolumeVariables& elemVolVars,
78 const SubControlVolumeFace& scvf,
79 const bool forceUpdateAll = false)
80 {
81 // fill the physics-related quantities of the caches
82 if (forceUpdateAll)
83 {
84 if constexpr (diffusionEnabled)
85 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
86 if constexpr (heatConductionEnabled)
87 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
88 }
89 else
90 {
91 if constexpr (diffusionEnabled && diffusionIsSolDependent)
92 fillDiffusion_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
93 if constexpr (heatConductionEnabled && heatConductionIsSolDependent)
94 fillHeatConduction_(scvfFluxVarsCache, element, fvGeometry, elemVolVars, scvf);
95 }
96 }
97
98private:
99
100 const Problem& problem() const
101 { return *problemPtr_; }
102
103
105 template<class FluxVariablesCache, class ElementVolumeVariables>
106 void fillDiffusion_(FluxVariablesCache& scvfFluxVarsCache,
107 const Element& element,
108 const FVElementGeometry& fvGeometry,
109 const ElementVolumeVariables& elemVolVars,
110 const SubControlVolumeFace& scvf)
111 {
112 using DiffusionType = typename ElementVolumeVariables::VolumeVariables::MolecularDiffusionType;
113 using DiffusionFiller = typename DiffusionType::Cache::Filler;
114 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
115
116 static constexpr int numPhases = ModelTraits::numFluidPhases();
117 static constexpr int numComponents = ModelTraits::numFluidComponents();
118
119 // forward to the filler of the diffusive quantities
120 if constexpr (FluidSystem::isTracerFluidSystem())
121 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
122 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
123 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *this);
124 else
125 for (unsigned int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
126 for (unsigned int compIdx = 0; compIdx < numComponents; ++compIdx)
127 if (compIdx != FluidSystem::getMainComponent(phaseIdx))
128 DiffusionFiller::fill(scvfFluxVarsCache, phaseIdx, compIdx, problem(), element, fvGeometry, elemVolVars, scvf, *this);
129 }
130
132 template<class FluxVariablesCache, class ElementVolumeVariables>
133 void fillHeatConduction_(FluxVariablesCache& scvfFluxVarsCache,
134 const Element& element,
135 const FVElementGeometry& fvGeometry,
136 const ElementVolumeVariables& elemVolVars,
137 const SubControlVolumeFace& scvf)
138 {
139 using HeatConductionType = typename ElementVolumeVariables::VolumeVariables::HeatConductionType;
140 using HeatConductionFiller = typename HeatConductionType::Cache::Filler;
141
142 // forward to the filler of the diffusive quantities
143 HeatConductionFiller::fill(scvfFluxVarsCache, problem(), element, fvGeometry, elemVolVars, scvf, *this);
144 }
145
146 const Problem* problemPtr_;
147};
148
149} // end namespace Dumux
150
151#endif
FreeFlowScalarFluxVariablesCacheFillerImplementation(const Problem &problem)
The constructor. Sets the problem pointer.
Definition: scalarfluxvariablescachefiller.hh:58
void fill(FluxVariablesCacheContainer &fluxVarsCacheContainer, FluxVariablesCache &scvfFluxVarsCache, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const bool forceUpdateAll=false)
function to fill the flux variables caches
Definition: scalarfluxvariablescachefiller.hh:73
Definition: scalarfluxvariablescachefiller.hh:27
Defines all properties used in Dumux.
Type traits for problem classes.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
std::decay_t< decltype(std::declval< Problem >().gridGeometry())> GridGeometry
Definition: common/typetraits/problem.hh:33