24#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_FOURIERS_LAW_HH
27#include <dune/common/dynvector.hh>
28#include <dune/common/dynmatrix.hh>
36template<
class TypeTag,
class DiscretizationMethod>
37class FouriersLawImplementation;
43template <
class TypeTag>
49 using Element =
typename GridView::template Codim<0>::Entity;
51 static constexpr int dim = GridView::dimension;
52 static constexpr int dimWorld = GridView::dimensionworld;
55 using FVElementGeometry =
typename GridGeometry::LocalView;
56 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
59 using ElementFluxVarsCache =
typename GridFluxVariablesCache::LocalView;
60 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
63 class MpfaFouriersLawCacheFiller
68 template<
class FluxVariablesCacheFiller>
69 static void fill(FluxVariablesCache& scvfFluxVarsCache,
70 const Problem& problem,
71 const Element& element,
72 const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars,
74 const SubControlVolumeFace& scvf,
75 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
78 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
79 scvfFluxVarsCache.updateHeatConduction(fluxVarsCacheFiller.secondaryInteractionVolume(),
80 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
81 fluxVarsCacheFiller.secondaryIvDataHandle());
83 scvfFluxVarsCache.updateHeatConduction(fluxVarsCacheFiller.primaryInteractionVolume(),
84 fluxVarsCacheFiller.primaryIvLocalFaceData(),
85 fluxVarsCacheFiller.primaryIvDataHandle());
90 class MpfaFouriersLawCache
93 using Stencil =
typename DualGridNodalIndexSet::NodalGridStencilType;
95 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
96 using PrimaryDataHandle =
typename ElementFluxVarsCache::PrimaryIvDataHandle::HeatConductionHandle;
97 using SecondaryDataHandle =
typename ElementFluxVarsCache::SecondaryIvDataHandle::HeatConductionHandle;
100 template<
bool doSecondary = cons
iderSecondaryIVs, std::enable_if_t<doSecondary,
int> = 0 >
101 void setHandlePointer_(
const SecondaryDataHandle& dataHandle)
102 { secondaryHandlePtr_ = &dataHandle; }
105 void setHandlePointer_(
const PrimaryDataHandle& dataHandle)
106 { primaryHandlePtr_ = &dataHandle; }
110 using Filler = MpfaFouriersLawCacheFiller;
120 template<
class IV,
class LocalFaceData,
class DataHandle>
121 void updateHeatConduction(
const IV& iv,
122 const LocalFaceData& localFaceData,
123 const DataHandle& dataHandle)
125 switchFluxSign_ = localFaceData.isOutsideFace();
126 stencil_ = &iv.stencil();
127 setHandlePointer_(dataHandle.heatConductionHandle());
131 const Stencil& heatConductionStencil()
const {
return *stencil_; }
134 const PrimaryDataHandle& heatConductionPrimaryDataHandle()
const {
return *primaryHandlePtr_; }
135 const SecondaryDataHandle& heatConductionSecondaryDataHandle()
const {
return *secondaryHandlePtr_; }
138 bool heatConductionSwitchFluxSign()
const {
return switchFluxSign_; }
141 bool switchFluxSign_;
144 const PrimaryDataHandle* primaryHandlePtr_;
145 const SecondaryDataHandle* secondaryHandlePtr_;
148 const Stencil* stencil_;
166 static Scalar
flux(
const Problem& problem,
167 const Element& element,
168 const FVElementGeometry& fvGeometry,
169 const ElementVolumeVariables& elemVolVars,
170 const SubControlVolumeFace& scvf,
171 const ElementFluxVarsCache& elemFluxVarsCache)
173 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
176 if (fluxVarsCache.usesSecondaryIv())
177 return flux_(problem, fluxVarsCache, fluxVarsCache.heatConductionSecondaryDataHandle());
179 return flux_(problem, fluxVarsCache, fluxVarsCache.heatConductionPrimaryDataHandle());
183 template<
class Problem,
class FluxVarsCache,
class DataHandle >
184 static Scalar flux_(
const Problem& problem,
185 const FluxVarsCache& cache,
186 const DataHandle& dataHandle)
188 const bool switchSign = cache.heatConductionSwitchFluxSign();
190 const auto localFaceIdx = cache.ivLocalFaceIndex();
191 const auto idxInOutside = cache.indexInOutsideFaces();
192 const auto& Tj = dataHandle.uj();
193 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
194 : (!switchSign ? dataHandle.T()[localFaceIdx]
195 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
196 Scalar scvfFlux = tij*Tj;
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:38
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Returns the heat flux within the porous medium (in J/s) across the given sub-control volume face.
Definition: flux/ccmpfa/fourierslaw.hh:166
MpfaFouriersLawCache Cache
Definition: flux/ccmpfa/fourierslaw.hh:157
Declares all properties used in Dumux.