12#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FOURIERS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_TPFA_FACET_COUPLING_FOURIERS_LAW_HH
18#include <dune/common/float_cmp.hh>
19#include <dune/common/fvector.hh>
34template<
class TypeTag,
44template<
class TypeTag>
53template<
class TypeTag>
62 using FVElementGeometry =
typename GridGeometry::LocalView;
63 using SubControlVolume =
typename GridGeometry::SubControlVolume;
64 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
67 using GridView =
typename GridGeometry::GridView;
68 using Element =
typename GridView::template Codim<0>::Entity;
73 class FacetCouplingFouriersLawCache
77 using Filler =
typename ParentType::Cache::Filler;
82 static constexpr int insideTijIdx = 0;
83 static constexpr int outsideTijIdx = 1;
84 static constexpr int facetTijIdx = 2;
87 using HeatConductionTransmissibilityContainer = std::array<Scalar, 3>;
90 template<
class Problem,
class ElementVolumeVariables >
91 void updateHeatConduction(
const Problem& problem,
92 const Element& element,
93 const FVElementGeometry& fvGeometry,
94 const ElementVolumeVariables& elemVolVars,
95 const SubControlVolumeFace &scvf)
97 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
104 Scalar heatConductionTij()
const
105 {
return tij_[insideTijIdx]; }
108 Scalar heatConductionTijInside()
const
109 {
return tij_[insideTijIdx]; }
112 Scalar heatConductionTijOutside()
const
113 {
return tij_[outsideTijIdx]; }
116 Scalar heatConductionTijFacet()
const
117 {
return tij_[facetTijIdx]; }
120 HeatConductionTransmissibilityContainer tij_;
125 using Cache = FacetCouplingFouriersLawCache;
132 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
133 static Scalar
flux(
const Problem& problem,
134 const Element& element,
135 const FVElementGeometry& fvGeometry,
136 const ElementVolumeVariables& elemVolVars,
137 const SubControlVolumeFace& scvf,
138 const ElementFluxVarsCache& elemFluxVarsCache)
140 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
141 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
144 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
145 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
146 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
147 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
150 const Scalar tInside = insideVolVars.temperature();
151 const Scalar tFacet = facetVolVars.temperature();
152 const Scalar tOutside = outsideVolVars.temperature();
154 Scalar flux = fluxVarsCache.heatConductionTijInside()*tInside
155 + fluxVarsCache.heatConductionTijFacet()*tFacet;
157 if (!scvf.boundary())
158 flux += fluxVarsCache.heatConductionTijOutside()*tOutside;
164 template<
class Problem,
class ElementVolumeVariables >
165 static typename Cache::HeatConductionTransmissibilityContainer
167 const Element& element,
168 const FVElementGeometry& fvGeometry,
169 const ElementVolumeVariables& elemVolVars,
170 const SubControlVolumeFace& scvf)
172 typename Cache::HeatConductionTransmissibilityContainer tij;
173 if (!problem.couplingManager().isCoupled(element, scvf))
176 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
181 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
183 const auto insideScvIdx = scvf.insideScvIdx();
184 const auto& insideScv = fvGeometry.scv(insideScvIdx);
185 const auto& insideVolVars = elemVolVars[insideScvIdx];
186 const auto wIn = Extrusion::area(fvGeometry, scvf)
188 insideVolVars.effectiveThermalConductivity(),
189 insideVolVars.extrusionFactor());
192 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
195 if (iBcTypes.hasOnlyNeumann())
197 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
198 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
199 /facetVolVars.extrusionFactor()
200 *
vtmv(scvf.unitOuterNormal(),
201 facetVolVars.effectiveThermalConductivity(),
202 scvf.unitOuterNormal());
211 if (!scvf.boundary())
213 const auto outsideScvIdx = scvf.outsideScvIdx();
214 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
215 const auto wOut = -1.0*Extrusion::area(fvGeometry, scvf)
217 outsideVolVars.effectiveThermalConductivity(),
218 outsideVolVars.extrusionFactor());
220 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
224 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
225 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
226 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
227 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
231 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
232 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
233 tij[Cache::outsideTijIdx] = 0.0;
238 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
239 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
240 tij[Cache::outsideTijIdx] = 0.0;
243 else if (iBcTypes.hasOnlyDirichlet())
245 tij[Cache::insideTijIdx] = wIn;
246 tij[Cache::outsideTijIdx] = 0.0;
247 tij[Cache::facetTijIdx] = -wIn;
250 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
260template<
class TypeTag>
269 using FVElementGeometry =
typename GridGeometry::LocalView;
270 using SubControlVolume =
typename GridGeometry::SubControlVolume;
271 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
274 using GridView =
typename GridGeometry::GridView;
275 using Element =
typename GridView::template Codim<0>::Entity;
280 class FacetCouplingFouriersLawCache
284 using Filler =
typename ParentType::Cache::Filler;
289 static constexpr int insideTijIdx = 0;
290 static constexpr int facetTijIdx = 1;
293 using HeatConductionTransmissibilityContainer = std::array<Scalar, 2>;
296 template<
class Problem,
class ElementVolumeVariables >
297 void updateHeatConduction(
const Problem& problem,
298 const Element& element,
299 const FVElementGeometry& fvGeometry,
300 const ElementVolumeVariables& elemVolVars,
301 const SubControlVolumeFace &scvf)
303 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
310 Scalar heatConductionTij()
const
311 {
return tij_[insideTijIdx]; }
314 Scalar heatConductionTijInside()
const
315 {
return tij_[insideTijIdx]; }
318 Scalar heatConductionTijFacet()
const
319 {
return tij_[facetTijIdx]; }
322 HeatConductionTransmissibilityContainer tij_;
327 using Cache = FacetCouplingFouriersLawCache;
334 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
335 static Scalar
flux(
const Problem& problem,
336 const Element& element,
337 const FVElementGeometry& fvGeometry,
338 const ElementVolumeVariables& elemVolVars,
339 const SubControlVolumeFace& scvf,
340 const ElementFluxVarsCache& elemFluxVarsCache)
342 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
343 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
346 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
347 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
348 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
350 return fluxVarsCache.heatConductionTijInside()*insideVolVars.temperature()
351 + fluxVarsCache.heatConductionTijFacet()*facetVolVars.temperature();
356 template<
class Problem,
class ElementVolumeVariables >
357 static typename Cache::HeatConductionTransmissibilityContainer
359 const Element& element,
360 const FVElementGeometry& fvGeometry,
361 const ElementVolumeVariables& elemVolVars,
362 const SubControlVolumeFace& scvf)
364 typename Cache::HeatConductionTransmissibilityContainer tij;
365 if (!problem.couplingManager().isCoupled(element, scvf))
368 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
373 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
378 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
379 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
381 const auto insideScvIdx = scvf.insideScvIdx();
382 const auto& insideScv = fvGeometry.scv(insideScvIdx);
383 const auto& insideVolVars = elemVolVars[insideScvIdx];
384 const auto wIn = Extrusion::area(fvGeometry, scvf)
386 insideVolVars.effectiveThermalConductivity(),
387 insideVolVars.extrusionFactor());
390 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
393 if (iBcTypes.hasOnlyNeumann())
398 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
399 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
400 /sqrt(facetVolVars.extrusionFactor())
401 *
vtmv(scvf.unitOuterNormal(),
402 facetVolVars.effectiveThermalConductivity(),
403 scvf.unitOuterNormal());
405 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
406 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
408 else if (iBcTypes.hasOnlyDirichlet())
410 tij[Cache::insideTijIdx] = wIn;
411 tij[Cache::facetTijIdx] = -wIn;
414 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:56
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:166
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:133
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:263
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:335
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:358
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:36
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fourierslaw.hh:34
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:100
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:26
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::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:36
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:880
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...