24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
37template<
class TypeTag,
class DiscretizationMethod>
38class FouriersLawImplementation;
44template <
class TypeTag>
51 using FVElementGeometry =
typename GridGeometry::LocalView;
52 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
54 using GridView =
typename GridGeometry::GridView;
56 using VolumeVariables =
typename ElementVolumeVariables::VolumeVariables;
57 using Element =
typename GridView::template Codim<0>::Entity;
59 using ElementFluxVarsCache =
typename GridFluxVariablesCache::LocalView;
60 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
62 static const int dim = GridView::dimension;
63 static const int dimWorld = GridView::dimensionworld;
66 class TpfaFouriersLawCacheFiller
71 template<
class FluxVariablesCacheFiller>
72 static void fill(FluxVariablesCache& scvfFluxVarsCache,
73 const Problem& problem,
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const SubControlVolumeFace& scvf,
78 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
80 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
85 class TpfaFouriersLawCache
88 using Filler = TpfaFouriersLawCacheFiller;
90 void updateHeatConduction(
const Problem& problem,
91 const Element& element,
92 const FVElementGeometry& fvGeometry,
93 const ElementVolumeVariables& elemVolVars,
94 const SubControlVolumeFace &scvf)
96 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
99 const Scalar& heatConductionTij()
const
125 static Scalar
flux(
const Problem& problem,
126 const Element& element,
127 const FVElementGeometry& fvGeometry,
128 const VolumeVariables& insideVolVars,
129 const VolumeVariables& outsideVolVars,
130 const SubControlVolumeFace& scvf,
131 const ElementFluxVarsCache& elemFluxVarsCache)
133 if constexpr (isMixedDimensional_)
134 if (scvf.numOutsideScv() != 1)
135 DUNE_THROW(Dune::Exception,
"This flux overload requires scvf.numOutsideScv() == 1");
138 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
141 const auto tInside = insideVolVars.temperature();
142 const auto tOutside = outsideVolVars.temperature();
144 return tij*(tInside - tOutside);
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)
162 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
165 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
166 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
167 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
169 return tij*(tInside - tOutside);
174 const Element& element,
175 const FVElementGeometry& fvGeometry,
176 const ElementVolumeVariables& elemVolVars,
177 const SubControlVolumeFace& scvf)
181 const auto insideScvIdx = scvf.insideScvIdx();
182 const auto& insideScv = fvGeometry.scv(insideScvIdx);
183 const auto& insideVolVars = elemVolVars[insideScvIdx];
185 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
189 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
191 tij = Extrusion::area(fvGeometry, scvf)*ti;
196 const auto outsideScvIdx = scvf.outsideScvIdx();
197 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
198 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
200 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
202 if constexpr (dim == dimWorld)
206 tj =
computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
212 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
221 static Scalar branchingFacetTemperature_(
const Problem& problem,
222 const Element& element,
223 const FVElementGeometry& fvGeometry,
224 const ElementVolumeVariables& elemVolVars,
225 const ElementFluxVarsCache& elemFluxVarsCache,
226 const SubControlVolumeFace& scvf,
227 Scalar insideTemperature,
230 Scalar sumTi(insideTi);
231 Scalar sumTempTi(insideTi*insideTemperature);
233 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
235 const auto outsideScvIdx = scvf.outsideScvIdx(i);
236 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
237 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
238 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
240 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
242 sumTempTi += outsideTi*outsideVolVars.temperature();
244 return sumTempTi/sumTi;
247 static constexpr bool isMixedDimensional_ =
static_cast<int>(GridView::dimension) <
static_cast<int>(GridView::dimensionworld);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename 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:48
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
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
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fourierslaw.hh:46
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:154
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:125
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:173
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:112
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...