12#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FOURIERS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_MPFA_FOURIERS_LAW_HH
15#include <dune/common/dynvector.hh>
16#include <dune/common/dynmatrix.hh>
24template<
class TypeTag,
class DiscretizationMethod>
25class FouriersLawImplementation;
31template <
class TypeTag>
37 using Element =
typename GridView::template Codim<0>::Entity;
39 static constexpr int dim = GridView::dimension;
40 static constexpr int dimWorld = GridView::dimensionworld;
43 using FVElementGeometry =
typename GridGeometry::LocalView;
44 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
47 using ElementFluxVarsCache =
typename GridFluxVariablesCache::LocalView;
48 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
51 class MpfaFouriersLawCacheFiller
56 template<
class FluxVariablesCacheFiller>
57 static void fill(FluxVariablesCache& scvfFluxVarsCache,
58 const Problem& problem,
59 const Element& element,
60 const FVElementGeometry& fvGeometry,
61 const ElementVolumeVariables& elemVolVars,
62 const SubControlVolumeFace& scvf,
63 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
66 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
67 scvfFluxVarsCache.updateHeatConduction(fluxVarsCacheFiller.secondaryInteractionVolume(),
68 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
69 fluxVarsCacheFiller.secondaryIvDataHandle());
71 scvfFluxVarsCache.updateHeatConduction(fluxVarsCacheFiller.primaryInteractionVolume(),
72 fluxVarsCacheFiller.primaryIvLocalFaceData(),
73 fluxVarsCacheFiller.primaryIvDataHandle());
78 class MpfaFouriersLawCache
81 using Stencil =
typename DualGridNodalIndexSet::NodalGridStencilType;
83 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
84 using PrimaryDataHandle =
typename ElementFluxVarsCache::PrimaryIvDataHandle::HeatConductionHandle;
85 using SecondaryDataHandle =
typename ElementFluxVarsCache::SecondaryIvDataHandle::HeatConductionHandle;
88 template<
bool doSecondary = cons
iderSecondaryIVs, std::enable_if_t<doSecondary,
int> = 0 >
89 void setHandlePointer_(
const SecondaryDataHandle& dataHandle)
90 { secondaryHandlePtr_ = &dataHandle; }
93 void setHandlePointer_(
const PrimaryDataHandle& dataHandle)
94 { primaryHandlePtr_ = &dataHandle; }
98 using Filler = MpfaFouriersLawCacheFiller;
108 template<
class IV,
class LocalFaceData,
class DataHandle>
109 void updateHeatConduction(
const IV& iv,
110 const LocalFaceData& localFaceData,
111 const DataHandle& dataHandle)
113 switchFluxSign_ = localFaceData.isOutsideFace();
114 stencil_ = &iv.stencil();
115 setHandlePointer_(dataHandle.heatConductionHandle());
119 const Stencil& heatConductionStencil()
const {
return *stencil_; }
122 const PrimaryDataHandle& heatConductionPrimaryDataHandle()
const {
return *primaryHandlePtr_; }
123 const SecondaryDataHandle& heatConductionSecondaryDataHandle()
const {
return *secondaryHandlePtr_; }
126 bool heatConductionSwitchFluxSign()
const {
return switchFluxSign_; }
129 bool switchFluxSign_;
132 const PrimaryDataHandle* primaryHandlePtr_;
133 const SecondaryDataHandle* secondaryHandlePtr_;
136 const Stencil* stencil_;
154 static Scalar
flux(
const Problem& problem,
155 const Element& element,
156 const FVElementGeometry& fvGeometry,
157 const ElementVolumeVariables& elemVolVars,
158 const SubControlVolumeFace& scvf,
159 const ElementFluxVarsCache& elemFluxVarsCache)
161 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
164 if (fluxVarsCache.usesSecondaryIv())
165 return flux_(problem, fluxVarsCache, fluxVarsCache.heatConductionSecondaryDataHandle());
167 return flux_(problem, fluxVarsCache, fluxVarsCache.heatConductionPrimaryDataHandle());
171 template<
class Problem,
class FluxVarsCache,
class DataHandle >
172 static Scalar flux_(
const Problem& problem,
173 const FluxVarsCache& cache,
174 const DataHandle& dataHandle)
176 const bool switchSign = cache.heatConductionSwitchFluxSign();
178 const auto localFaceIdx = cache.ivLocalFaceIndex();
179 const auto idxInOutside = cache.indexInOutsideFaces();
180 const auto& Tj = dataHandle.uj();
181 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
182 : (!switchSign ? dataHandle.T()[localFaceIdx]
183 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
184 Scalar scvfFlux = tij*Tj;
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:154
MpfaFouriersLawCache Cache
Definition: flux/ccmpfa/fourierslaw.hh:145
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:26
Defines all properties used in Dumux.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.