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;
61 static const int dim = GridView::dimension;
62 static const int dimWorld = GridView::dimensionworld;
65 class TpfaFouriersLawCacheFiller
70 template<
class FluxVariablesCacheFiller>
71 static void fill(FluxVariablesCache& scvfFluxVarsCache,
72 const Problem& problem,
73 const Element& element,
74 const FVElementGeometry& fvGeometry,
75 const ElementVolumeVariables& elemVolVars,
76 const SubControlVolumeFace& scvf,
77 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
79 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
84 class TpfaFouriersLawCache
87 using Filler = TpfaFouriersLawCacheFiller;
89 void updateHeatConduction(
const Problem& problem,
90 const Element& element,
91 const FVElementGeometry& fvGeometry,
92 const ElementVolumeVariables& elemVolVars,
93 const SubControlVolumeFace &scvf)
95 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
98 const Scalar& heatConductionTij()
const
124 static Scalar
flux(
const Problem& problem,
125 const Element& element,
126 const FVElementGeometry& fvGeometry,
127 const VolumeVariables& insideVolVars,
128 const VolumeVariables& outsideVolVars,
129 const SubControlVolumeFace& scvf,
130 const ElementFluxVarsCache& elemFluxVarsCache)
132 if constexpr (isMixedDimensional_)
133 if (scvf.numOutsideScv() != 1)
134 DUNE_THROW(Dune::Exception,
"This flux overload requires scvf.numOutsideScv() == 1");
137 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
140 const auto tInside = insideVolVars.temperature();
141 const auto tOutside = outsideVolVars.temperature();
143 return tij*(tInside - tOutside);
153 static Scalar
flux(
const Problem& problem,
154 const Element& element,
155 const FVElementGeometry& fvGeometry,
156 const ElementVolumeVariables& elemVolVars,
157 const SubControlVolumeFace& scvf,
158 const ElementFluxVarsCache& elemFluxVarsCache)
161 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
164 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
165 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
166 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
168 return tij*(tInside - tOutside);
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const SubControlVolumeFace& scvf)
180 const auto insideScvIdx = scvf.insideScvIdx();
181 const auto& insideScv = fvGeometry.scv(insideScvIdx);
182 const auto& insideVolVars = elemVolVars[insideScvIdx];
184 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
188 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
190 tij = Extrusion::area(scvf)*ti;
195 const auto outsideScvIdx = scvf.outsideScvIdx();
196 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
197 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
199 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
201 if constexpr (dim == dimWorld)
211 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
220 static Scalar branchingFacetTemperature_(
const Problem& problem,
221 const Element& element,
222 const FVElementGeometry& fvGeometry,
223 const ElementVolumeVariables& elemVolVars,
224 const ElementFluxVarsCache& elemFluxVarsCache,
225 const SubControlVolumeFace& scvf,
226 Scalar insideTemperature,
229 Scalar sumTi(insideTi);
230 Scalar sumTempTi(insideTi*insideTemperature);
232 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
234 const auto outsideScvIdx = scvf.outsideScvIdx(i);
235 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
236 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
237 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
239 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
241 sumTempTi += outsideTi*outsideVolVars.temperature();
243 return sumTempTi/sumTi;
246 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:47
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
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:153
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:124
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:172
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:111
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...