24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
37template<
class TypeTag, DiscretizationMethod discMethod>
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 Element =
typename GridView::template Codim<0>::Entity;
60 static const int dim = GridView::dimension;
61 static const int dimWorld = GridView::dimensionworld;
64 class TpfaFouriersLawCacheFiller
69 template<
class FluxVariablesCacheFiller>
70 static void fill(FluxVariablesCache& scvfFluxVarsCache,
71 const Problem& problem,
72 const Element& element,
73 const FVElementGeometry& fvGeometry,
74 const ElementVolumeVariables& elemVolVars,
75 const SubControlVolumeFace& scvf,
76 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
78 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
83 class TpfaFouriersLawCache
86 using Filler = TpfaFouriersLawCacheFiller;
88 void updateHeatConduction(
const Problem& problem,
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const SubControlVolumeFace &scvf)
94 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
97 const Scalar& heatConductionTij()
const
118 static Scalar
flux(
const Problem& problem,
119 const Element& element,
120 const FVElementGeometry& fvGeometry,
121 const ElementVolumeVariables& elemVolVars,
122 const SubControlVolumeFace& scvf,
123 const ElementFluxVarsCache& elemFluxVarsCache)
126 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
129 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
130 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
131 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
133 return tij*(tInside - tOutside);
138 const Element& element,
139 const FVElementGeometry& fvGeometry,
140 const ElementVolumeVariables& elemVolVars,
141 const SubControlVolumeFace& scvf)
145 const auto insideScvIdx = scvf.insideScvIdx();
146 const auto& insideScv = fvGeometry.scv(insideScvIdx);
147 const auto& insideVolVars = elemVolVars[insideScvIdx];
149 const auto insideLambda = insideVolVars.effectiveThermalConductivity();
153 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
155 tij = Extrusion::area(scvf)*ti;
160 const auto outsideScvIdx = scvf.outsideScvIdx();
161 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
162 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
164 const auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
176 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
185 static Scalar branchingFacetTemperature_(
const Problem& problem,
186 const Element& element,
187 const FVElementGeometry& fvGeometry,
188 const ElementVolumeVariables& elemVolVars,
189 const ElementFluxVarsCache& elemFluxVarsCache,
190 const SubControlVolumeFace& scvf,
191 Scalar insideTemperature,
194 Scalar sumTi(insideTi);
195 Scalar sumTempTi(insideTi*insideTemperature);
197 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
199 const auto outsideScvIdx = scvf.outsideScvIdx(i);
200 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
201 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
202 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
204 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
206 sumTempTi += outsideTi*outsideVolVars.temperature();
208 return sumTempTi/sumTi;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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 (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declaration of the method-specific implementation
Definition: flux/fourierslaw.hh:37
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fourierslaw.hh:46
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:137
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:118
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:109
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...