24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FOURIERS_LAW_HH
30#include <dune/common/float_cmp.hh>
31#include <dune/common/fvector.hh>
45template<
class TypeTag,
55template<
class TypeTag>
64template<
class TypeTag>
73 using FVElementGeometry =
typename GridGeometry::LocalView;
74 using SubControlVolume =
typename GridGeometry::SubControlVolume;
75 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
77 using GridView =
typename GridGeometry::GridView;
78 using Element =
typename GridView::template Codim<0>::Entity;
83 class FacetCouplingFouriersLawCache
87 using Filler =
typename ParentType::Cache::Filler;
92 static constexpr int insideTijIdx = 0;
93 static constexpr int outsideTijIdx = 1;
94 static constexpr int facetTijIdx = 2;
97 using HeatConductionTransmissibilityContainer = std::array<Scalar, 3>;
100 template<
class Problem,
class ElementVolumeVariables >
101 void updateHeatConduction(
const Problem& problem,
102 const Element& element,
103 const FVElementGeometry& fvGeometry,
104 const ElementVolumeVariables& elemVolVars,
105 const SubControlVolumeFace &scvf)
107 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
114 Scalar heatConductionTij()
const
115 {
return tij_[insideTijIdx]; }
118 Scalar heatConductionTijInside()
const
119 {
return tij_[insideTijIdx]; }
122 Scalar heatConductionTijOutside()
const
123 {
return tij_[outsideTijIdx]; }
126 Scalar heatConductionTijFacet()
const
127 {
return tij_[facetTijIdx]; }
130 HeatConductionTransmissibilityContainer tij_;
135 using Cache = FacetCouplingFouriersLawCache;
141 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
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)
149 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
150 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
153 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
154 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
155 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
156 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
159 const Scalar tInside = insideVolVars.temperature();
160 const Scalar tFacet = facetVolVars.temperature();
161 const Scalar tOutside = outsideVolVars.temperature();
163 Scalar flux = fluxVarsCache.heatConductionTijInside()*tInside
164 + fluxVarsCache.heatConductionTijFacet()*tFacet;
166 if (!scvf.boundary())
167 flux += fluxVarsCache.heatConductionTijOutside()*tOutside;
173 template<
class Problem,
class ElementVolumeVariables >
174 static typename Cache::HeatConductionTransmissibilityContainer
176 const Element& element,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const SubControlVolumeFace& scvf)
181 typename Cache::HeatConductionTransmissibilityContainer tij;
182 if (!problem.couplingManager().isCoupled(element, scvf))
185 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
190 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
192 const auto insideScvIdx = scvf.insideScvIdx();
193 const auto& insideScv = fvGeometry.scv(insideScvIdx);
194 const auto& insideVolVars = elemVolVars[insideScvIdx];
197 insideVolVars.effectiveThermalConductivity(),
198 insideVolVars.extrusionFactor());
201 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
204 if (iBcTypes.hasOnlyNeumann())
206 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
207 const auto wFacet = 2.0*scvf.area()*insideVolVars.extrusionFactor()
208 /facetVolVars.extrusionFactor()
209 *
vtmv(scvf.unitOuterNormal(),
210 facetVolVars.effectiveThermalConductivity(),
211 scvf.unitOuterNormal());
220 if (!scvf.boundary())
222 const auto outsideScvIdx = scvf.outsideScvIdx();
223 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
225 fvGeometry.scv(outsideScvIdx),
226 outsideVolVars.effectiveThermalConductivity(),
227 outsideVolVars.extrusionFactor());
229 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
233 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
234 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
235 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
236 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
240 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
241 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
242 tij[Cache::outsideTijIdx] = 0.0;
247 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
248 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
249 tij[Cache::outsideTijIdx] = 0.0;
252 else if (iBcTypes.hasOnlyDirichlet())
254 tij[Cache::insideTijIdx] = wIn;
255 tij[Cache::outsideTijIdx] = 0.0;
256 tij[Cache::facetTijIdx] = -wIn;
259 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
269template<
class TypeTag>
278 using FVElementGeometry =
typename GridGeometry::LocalView;
279 using SubControlVolume =
typename GridGeometry::SubControlVolume;
280 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
282 using GridView =
typename GridGeometry::GridView;
283 using Element =
typename GridView::template Codim<0>::Entity;
288 class FacetCouplingFouriersLawCache
292 using Filler =
typename ParentType::Cache::Filler;
297 static constexpr int insideTijIdx = 0;
298 static constexpr int facetTijIdx = 1;
301 using HeatConductionTransmissibilityContainer = std::array<Scalar, 2>;
304 template<
class Problem,
class ElementVolumeVariables >
305 void updateHeatConduction(
const Problem& problem,
306 const Element& element,
307 const FVElementGeometry& fvGeometry,
308 const ElementVolumeVariables& elemVolVars,
309 const SubControlVolumeFace &scvf)
311 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
318 Scalar heatConductionTij()
const
319 {
return tij_[insideTijIdx]; }
322 Scalar heatConductionTijInside()
const
323 {
return tij_[insideTijIdx]; }
326 Scalar heatConductionTijFacet()
const
327 {
return tij_[facetTijIdx]; }
330 HeatConductionTransmissibilityContainer tij_;
335 using Cache = FacetCouplingFouriersLawCache;
341 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
342 static Scalar
flux(
const Problem& problem,
343 const Element& element,
344 const FVElementGeometry& fvGeometry,
345 const ElementVolumeVariables& elemVolVars,
346 const SubControlVolumeFace& scvf,
347 const ElementFluxVarsCache& elemFluxVarsCache)
349 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
350 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
353 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
354 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
355 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
357 return fluxVarsCache.heatConductionTijInside()*insideVolVars.temperature()
358 + fluxVarsCache.heatConductionTijFacet()*facetVolVars.temperature();
363 template<
class Problem,
class ElementVolumeVariables >
364 static typename Cache::HeatConductionTransmissibilityContainer
366 const Element& element,
367 const FVElementGeometry& fvGeometry,
368 const ElementVolumeVariables& elemVolVars,
369 const SubControlVolumeFace& scvf)
371 typename Cache::HeatConductionTransmissibilityContainer tij;
372 if (!problem.couplingManager().isCoupled(element, scvf))
375 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
380 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
385 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
386 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
388 const auto insideScvIdx = scvf.insideScvIdx();
389 const auto& insideScv = fvGeometry.scv(insideScvIdx);
390 const auto& insideVolVars = elemVolVars[insideScvIdx];
393 insideVolVars.effectiveThermalConductivity(),
394 insideVolVars.extrusionFactor());
397 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
400 if (iBcTypes.hasOnlyNeumann())
405 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
406 const auto wFacet = 2.0*scvf.area()*insideVolVars.extrusionFactor()
407 /sqrt(facetVolVars.extrusionFactor())
408 *
vtmv(scvf.unitOuterNormal(),
409 facetVolVars.effectiveThermalConductivity(),
410 scvf.unitOuterNormal());
412 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
413 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
415 else if (iBcTypes.hasOnlyDirichlet())
417 tij[Cache::insideTijIdx] = wIn;
418 tij[Cache::facetTijIdx] = -wIn;
421 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
Define some often used mathematical functions.
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
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:840
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:47
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:108
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:47
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:67
static Cache::HeatConductionTransmissibilityContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:175
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the conductive heat flux.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:142
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:272
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Compute the conductive heat flux.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:342
static Cache::HeatConductionTransmissibilityContainer calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:365
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.