12#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
25template<
class TypeTag,
class DiscretizationMethod>
26class FouriersLawImplementation;
32template <
class TypeTag>
39 using FVElementGeometry =
typename GridGeometry::LocalView;
40 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
42 using GridView =
typename GridGeometry::GridView;
44 using VolumeVariables =
typename ElementVolumeVariables::VolumeVariables;
45 using Element =
typename GridView::template Codim<0>::Entity;
47 using ElementFluxVarsCache =
typename GridFluxVariablesCache::LocalView;
48 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
50 static const int dim = GridView::dimension;
51 static const int dimWorld = GridView::dimensionworld;
54 class TpfaFouriersLawCacheFiller
59 template<
class FluxVariablesCacheFiller>
60 static void fill(FluxVariablesCache& scvfFluxVarsCache,
61 const Problem& problem,
62 const Element& element,
63 const FVElementGeometry& fvGeometry,
64 const ElementVolumeVariables& elemVolVars,
65 const SubControlVolumeFace& scvf,
66 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
68 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
73 class TpfaFouriersLawCache
76 using Filler = TpfaFouriersLawCacheFiller;
78 void updateHeatConduction(
const Problem& problem,
79 const Element& element,
80 const FVElementGeometry& fvGeometry,
81 const ElementVolumeVariables& elemVolVars,
82 const SubControlVolumeFace &scvf)
84 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
87 const Scalar& heatConductionTij()
const
113 static Scalar
flux(
const Problem& problem,
114 const Element& element,
115 const FVElementGeometry& fvGeometry,
116 const VolumeVariables& insideVolVars,
117 const VolumeVariables& outsideVolVars,
118 const SubControlVolumeFace& scvf,
119 const ElementFluxVarsCache& elemFluxVarsCache)
121 if constexpr (isMixedDimensional_)
122 if (scvf.numOutsideScv() != 1)
123 DUNE_THROW(Dune::Exception,
"This flux overload requires scvf.numOutsideScv() == 1");
126 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
129 const auto tInside = insideVolVars.temperature();
130 const auto tOutside = outsideVolVars.temperature();
132 return tij*(tInside - tOutside);
142 static Scalar
flux(
const Problem& problem,
143 const Element& element,
144 const FVElementGeometry& fvGeometry,
145 const ElementVolumeVariables& elemVolVars,
146 const SubControlVolumeFace& scvf,
147 const ElementFluxVarsCache& elemFluxVarsCache)
150 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
153 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
154 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
155 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
157 return tij*(tInside - tOutside);
162 const Element& element,
163 const FVElementGeometry& fvGeometry,
164 const ElementVolumeVariables& elemVolVars,
165 const SubControlVolumeFace& scvf)
169 const auto insideScvIdx = scvf.insideScvIdx();
170 const auto& insideScv = fvGeometry.scv(insideScvIdx);
171 const auto& insideVolVars = elemVolVars[insideScvIdx];
173 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
177 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
179 tij = Extrusion::area(fvGeometry, scvf)*ti;
184 const auto outsideScvIdx = scvf.outsideScvIdx();
185 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
186 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
188 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
190 if constexpr (dim == dimWorld)
194 tj =
computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
200 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
209 static Scalar branchingFacetTemperature_(
const Problem& problem,
210 const Element& element,
211 const FVElementGeometry& fvGeometry,
212 const ElementVolumeVariables& elemVolVars,
213 const ElementFluxVarsCache& elemFluxVarsCache,
214 const SubControlVolumeFace& scvf,
215 Scalar insideTemperature,
218 Scalar sumTi(insideTi);
219 Scalar sumTempTi(insideTi*insideTemperature);
221 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
223 const auto outsideScvIdx = scvf.outsideScvIdx(i);
224 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
225 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
226 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
228 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
230 sumTempTi += outsideTi*outsideVolVars.temperature();
232 return sumTempTi/sumTi;
235 static constexpr bool isMixedDimensional_ =
static_cast<int>(GridView::dimension) <
static_cast<int>(GridView::dimensionworld);
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fourierslaw.hh:34
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/cctpfa/fourierslaw.hh:142
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &insideVolVars, const VolumeVariables &outsideVolVars, 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/cctpfa/fourierslaw.hh:113
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Compute transmissibilities.
Definition: flux/cctpfa/fourierslaw.hh:161
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:100
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:26
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::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:36
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...