24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_HH
38template<
class TypeTag, DiscretizationMethod discMethod>
39class FouriersLawImplementation;
45template <
class TypeTag>
52 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
53 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
56 using Element =
typename GridView::template Codim<0>::Entity;
60 static const int dim = GridView::dimension;
61 static const int dimWorld = GridView::dimensionworld;
63 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
68 class TpfaFouriersLawCacheFiller
73 template<
class FluxVariablesCacheFiller>
74 static void fill(FluxVariablesCache& scvfFluxVarsCache,
75 const Problem& problem,
76 const Element& element,
77 const FVElementGeometry& fvGeometry,
78 const ElementVolumeVariables& elemVolVars,
79 const SubControlVolumeFace& scvf,
80 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
82 scvfFluxVarsCache.updateHeatConduction(problem, element, fvGeometry, elemVolVars, scvf);
87 class TpfaFouriersLawCache
90 using Filler = TpfaFouriersLawCacheFiller;
92 void updateHeatConduction(
const Problem& problem,
93 const Element& element,
94 const FVElementGeometry& fvGeometry,
95 const ElementVolumeVariables& elemVolVars,
96 const SubControlVolumeFace &scvf)
98 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
101 const Scalar& heatConductionTij()
const
116 static Scalar
flux(
const Problem& problem,
117 const Element& element,
118 const FVElementGeometry& fvGeometry,
119 const ElementVolumeVariables& elemVolVars,
120 const SubControlVolumeFace& scvf,
121 const ElementFluxVarsCache& elemFluxVarsCache)
124 Scalar tij = elemFluxVarsCache[scvf].heatConductionTij();
127 const auto tInside = elemVolVars[scvf.insideScvIdx()].temperature();
128 const auto tOutside = scvf.numOutsideScvs() == 1 ? elemVolVars[scvf.outsideScvIdx()].temperature()
129 : branchingFacetTemperature_(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf, tInside, tij);
131 return tij*(tInside - tOutside);
136 const Element& element,
137 const FVElementGeometry& fvGeometry,
138 const ElementVolumeVariables& elemVolVars,
139 const SubControlVolumeFace& scvf)
143 const auto insideScvIdx = scvf.insideScvIdx();
144 const auto& insideScv = fvGeometry.scv(insideScvIdx);
145 const auto& insideVolVars = elemVolVars[insideScvIdx];
147 const auto insideLambda = Deprecated::template effectiveThermalConductivity<ThermalConductivityModel>(
148 insideVolVars, problem.spatialParams(), element, fvGeometry, insideScv);
152 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
154 tij = scvf.area()*ti;
159 const auto outsideScvIdx = scvf.outsideScvIdx();
160 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
161 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
162 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
164 const auto outsideLambda = Deprecated::template effectiveThermalConductivity<ThermalConductivityModel>(
165 outsideVolVars, problem.spatialParams(), outsideElement, fvGeometry, outsideScv);
177 tij = scvf.area()*(ti * tj)/(ti + tj);
186 static Scalar branchingFacetTemperature_(
const Problem& problem,
187 const Element& element,
188 const FVElementGeometry& fvGeometry,
189 const ElementVolumeVariables& elemVolVars,
190 const ElementFluxVarsCache& elemFluxVarsCache,
191 const SubControlVolumeFace& scvf,
192 Scalar insideTemperature,
195 Scalar sumTi(insideTi);
196 Scalar sumTempTi(insideTi*insideTemperature);
198 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
200 const auto outsideScvIdx = scvf.outsideScvIdx(i);
201 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
202 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
203 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
205 auto outsideTi = outsideFluxVarsCache.heatConductionTij();
207 sumTempTi += outsideTi*outsideVolVars.temperature();
209 return sumTempTi/sumTi;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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: fourierslaw.hh:37
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: cctpfa/fourierslaw.hh:47
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Compute transmissibilities.
Definition: cctpfa/fourierslaw.hh:135
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the heat condution flux assuming thermal equilibrium.
Definition: cctpfa/fourierslaw.hh:116
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: cctpfa/fourierslaw.hh:113
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...