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>
46template<
class TypeTag,
56template<
class TypeTag>
65template<
class TypeTag>
74 using FVElementGeometry =
typename GridGeometry::LocalView;
75 using SubControlVolume =
typename GridGeometry::SubControlVolume;
76 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
79 using GridView =
typename GridGeometry::GridView;
80 using Element =
typename GridView::template Codim<0>::Entity;
85 class FacetCouplingFouriersLawCache
89 using Filler =
typename ParentType::Cache::Filler;
94 static constexpr int insideTijIdx = 0;
95 static constexpr int outsideTijIdx = 1;
96 static constexpr int facetTijIdx = 2;
99 using HeatConductionTransmissibilityContainer = std::array<Scalar, 3>;
102 template<
class Problem,
class ElementVolumeVariables >
103 void updateHeatConduction(
const Problem& problem,
104 const Element& element,
105 const FVElementGeometry& fvGeometry,
106 const ElementVolumeVariables& elemVolVars,
107 const SubControlVolumeFace &scvf)
109 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
116 Scalar heatConductionTij()
const
117 {
return tij_[insideTijIdx]; }
120 Scalar heatConductionTijInside()
const
121 {
return tij_[insideTijIdx]; }
124 Scalar heatConductionTijOutside()
const
125 {
return tij_[outsideTijIdx]; }
128 Scalar heatConductionTijFacet()
const
129 {
return tij_[facetTijIdx]; }
132 HeatConductionTransmissibilityContainer tij_;
137 using Cache = FacetCouplingFouriersLawCache;
144 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
145 static Scalar
flux(
const Problem& problem,
146 const Element& element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& elemVolVars,
149 const SubControlVolumeFace& scvf,
150 const ElementFluxVarsCache& elemFluxVarsCache)
152 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
153 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
156 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
157 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
158 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
159 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
162 const Scalar tInside = insideVolVars.temperature();
163 const Scalar tFacet = facetVolVars.temperature();
164 const Scalar tOutside = outsideVolVars.temperature();
166 Scalar flux = fluxVarsCache.heatConductionTijInside()*tInside
167 + fluxVarsCache.heatConductionTijFacet()*tFacet;
169 if (!scvf.boundary())
170 flux += fluxVarsCache.heatConductionTijOutside()*tOutside;
176 template<
class Problem,
class ElementVolumeVariables >
177 static typename Cache::HeatConductionTransmissibilityContainer
179 const Element& element,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars,
182 const SubControlVolumeFace& scvf)
184 typename Cache::HeatConductionTransmissibilityContainer tij;
185 if (!problem.couplingManager().isCoupled(element, scvf))
188 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
193 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
195 const auto insideScvIdx = scvf.insideScvIdx();
196 const auto& insideScv = fvGeometry.scv(insideScvIdx);
197 const auto& insideVolVars = elemVolVars[insideScvIdx];
198 const auto wIn = Extrusion::area(fvGeometry, scvf)
200 insideVolVars.effectiveThermalConductivity(),
201 insideVolVars.extrusionFactor());
204 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
207 if (iBcTypes.hasOnlyNeumann())
209 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
210 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
211 /facetVolVars.extrusionFactor()
212 *
vtmv(scvf.unitOuterNormal(),
213 facetVolVars.effectiveThermalConductivity(),
214 scvf.unitOuterNormal());
223 if (!scvf.boundary())
225 const auto outsideScvIdx = scvf.outsideScvIdx();
226 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
227 const auto wOut = -1.0*Extrusion::area(fvGeometry, scvf)
229 outsideVolVars.effectiveThermalConductivity(),
230 outsideVolVars.extrusionFactor());
232 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
236 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
237 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
238 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
239 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
243 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
244 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
245 tij[Cache::outsideTijIdx] = 0.0;
250 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
251 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
252 tij[Cache::outsideTijIdx] = 0.0;
255 else if (iBcTypes.hasOnlyDirichlet())
257 tij[Cache::insideTijIdx] = wIn;
258 tij[Cache::outsideTijIdx] = 0.0;
259 tij[Cache::facetTijIdx] = -wIn;
262 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
272template<
class TypeTag>
281 using FVElementGeometry =
typename GridGeometry::LocalView;
282 using SubControlVolume =
typename GridGeometry::SubControlVolume;
283 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
286 using GridView =
typename GridGeometry::GridView;
287 using Element =
typename GridView::template Codim<0>::Entity;
292 class FacetCouplingFouriersLawCache
296 using Filler =
typename ParentType::Cache::Filler;
301 static constexpr int insideTijIdx = 0;
302 static constexpr int facetTijIdx = 1;
305 using HeatConductionTransmissibilityContainer = std::array<Scalar, 2>;
308 template<
class Problem,
class ElementVolumeVariables >
309 void updateHeatConduction(
const Problem& problem,
310 const Element& element,
311 const FVElementGeometry& fvGeometry,
312 const ElementVolumeVariables& elemVolVars,
313 const SubControlVolumeFace &scvf)
315 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
322 Scalar heatConductionTij()
const
323 {
return tij_[insideTijIdx]; }
326 Scalar heatConductionTijInside()
const
327 {
return tij_[insideTijIdx]; }
330 Scalar heatConductionTijFacet()
const
331 {
return tij_[facetTijIdx]; }
334 HeatConductionTransmissibilityContainer tij_;
339 using Cache = FacetCouplingFouriersLawCache;
346 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
347 static Scalar
flux(
const Problem& problem,
348 const Element& element,
349 const FVElementGeometry& fvGeometry,
350 const ElementVolumeVariables& elemVolVars,
351 const SubControlVolumeFace& scvf,
352 const ElementFluxVarsCache& elemFluxVarsCache)
354 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
355 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
358 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
359 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
360 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
362 return fluxVarsCache.heatConductionTijInside()*insideVolVars.temperature()
363 + fluxVarsCache.heatConductionTijFacet()*facetVolVars.temperature();
368 template<
class Problem,
class ElementVolumeVariables >
369 static typename Cache::HeatConductionTransmissibilityContainer
371 const Element& element,
372 const FVElementGeometry& fvGeometry,
373 const ElementVolumeVariables& elemVolVars,
374 const SubControlVolumeFace& scvf)
376 typename Cache::HeatConductionTransmissibilityContainer tij;
377 if (!problem.couplingManager().isCoupled(element, scvf))
380 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
385 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
390 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
391 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
393 const auto insideScvIdx = scvf.insideScvIdx();
394 const auto& insideScv = fvGeometry.scv(insideScvIdx);
395 const auto& insideVolVars = elemVolVars[insideScvIdx];
396 const auto wIn = Extrusion::area(fvGeometry, scvf)
398 insideVolVars.effectiveThermalConductivity(),
399 insideVolVars.extrusionFactor());
402 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
405 if (iBcTypes.hasOnlyNeumann())
410 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
411 const auto wFacet = 2.0*Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor()
412 /sqrt(facetVolVars.extrusionFactor())
413 *
vtmv(scvf.unitOuterNormal(),
414 facetVolVars.effectiveThermalConductivity(),
415 scvf.unitOuterNormal());
417 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
418 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
420 else if (iBcTypes.hasOnlyDirichlet())
422 tij[Cache::insideTijIdx] = wIn;
423 tij[Cache::facetTijIdx] = -wIn;
426 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
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:48
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:863
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
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
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:112
Forward declaration of the implementation.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:48
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim=dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:68
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:178
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:145
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:275
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:347
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:370
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.