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;
143 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
144 static Scalar
flux(
const Problem& problem,
145 const Element& element,
146 const FVElementGeometry& fvGeometry,
147 const ElementVolumeVariables& elemVolVars,
148 const SubControlVolumeFace& scvf,
149 const ElementFluxVarsCache& elemFluxVarsCache)
151 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
152 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
155 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
156 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
157 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
158 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
161 const Scalar tInside = insideVolVars.temperature();
162 const Scalar tFacet = facetVolVars.temperature();
163 const Scalar tOutside = outsideVolVars.temperature();
165 Scalar flux = fluxVarsCache.heatConductionTijInside()*tInside
166 + fluxVarsCache.heatConductionTijFacet()*tFacet;
168 if (!scvf.boundary())
169 flux += fluxVarsCache.heatConductionTijOutside()*tOutside;
175 template<
class Problem,
class ElementVolumeVariables >
176 static typename Cache::HeatConductionTransmissibilityContainer
178 const Element& element,
179 const FVElementGeometry& fvGeometry,
180 const ElementVolumeVariables& elemVolVars,
181 const SubControlVolumeFace& scvf)
183 typename Cache::HeatConductionTransmissibilityContainer tij;
184 if (!problem.couplingManager().isCoupled(element, scvf))
187 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
192 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
194 const auto insideScvIdx = scvf.insideScvIdx();
195 const auto& insideScv = fvGeometry.scv(insideScvIdx);
196 const auto& insideVolVars = elemVolVars[insideScvIdx];
197 const auto wIn = Extrusion::area(scvf)
199 insideVolVars.effectiveThermalConductivity(),
200 insideVolVars.extrusionFactor());
203 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
206 if (iBcTypes.hasOnlyNeumann())
208 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
209 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
210 /facetVolVars.extrusionFactor()
211 *
vtmv(scvf.unitOuterNormal(),
212 facetVolVars.effectiveThermalConductivity(),
213 scvf.unitOuterNormal());
222 if (!scvf.boundary())
224 const auto outsideScvIdx = scvf.outsideScvIdx();
225 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
226 const auto wOut = -1.0*Extrusion::area(scvf)
228 outsideVolVars.effectiveThermalConductivity(),
229 outsideVolVars.extrusionFactor());
231 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
235 const Scalar factor = wIn * wFacet / ( wIn * wOut * ( 2.0 * xi - 1.0 ) + wFacet * ( xi * ( wIn + wOut ) + wFacet ) );
236 tij[Cache::insideTijIdx] = factor * ( wOut * xi + wFacet );
237 tij[Cache::outsideTijIdx] = factor * ( wOut * ( 1.0 - xi ) );
238 tij[Cache::facetTijIdx] = factor * ( - wOut - wFacet );
242 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
243 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
244 tij[Cache::outsideTijIdx] = 0.0;
249 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
250 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
251 tij[Cache::outsideTijIdx] = 0.0;
254 else if (iBcTypes.hasOnlyDirichlet())
256 tij[Cache::insideTijIdx] = wIn;
257 tij[Cache::outsideTijIdx] = 0.0;
258 tij[Cache::facetTijIdx] = -wIn;
261 DUNE_THROW(Dune::NotImplemented,
"Interior boundary types other than pure Dirichlet or Neumann");
271template<
class TypeTag>
280 using FVElementGeometry =
typename GridGeometry::LocalView;
281 using SubControlVolume =
typename GridGeometry::SubControlVolume;
282 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
285 using GridView =
typename GridGeometry::GridView;
286 using Element =
typename GridView::template Codim<0>::Entity;
291 class FacetCouplingFouriersLawCache
295 using Filler =
typename ParentType::Cache::Filler;
300 static constexpr int insideTijIdx = 0;
301 static constexpr int facetTijIdx = 1;
304 using HeatConductionTransmissibilityContainer = std::array<Scalar, 2>;
307 template<
class Problem,
class ElementVolumeVariables >
308 void updateHeatConduction(
const Problem& problem,
309 const Element& element,
310 const FVElementGeometry& fvGeometry,
311 const ElementVolumeVariables& elemVolVars,
312 const SubControlVolumeFace &scvf)
314 tij_ = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
321 Scalar heatConductionTij()
const
322 {
return tij_[insideTijIdx]; }
325 Scalar heatConductionTijInside()
const
326 {
return tij_[insideTijIdx]; }
329 Scalar heatConductionTijFacet()
const
330 {
return tij_[facetTijIdx]; }
333 HeatConductionTransmissibilityContainer tij_;
338 using Cache = FacetCouplingFouriersLawCache;
344 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache >
345 static Scalar
flux(
const Problem& problem,
346 const Element& element,
347 const FVElementGeometry& fvGeometry,
348 const ElementVolumeVariables& elemVolVars,
349 const SubControlVolumeFace& scvf,
350 const ElementFluxVarsCache& elemFluxVarsCache)
352 if (!problem.couplingManager().isOnInteriorBoundary(element, scvf))
353 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
356 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
357 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
358 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
360 return fluxVarsCache.heatConductionTijInside()*insideVolVars.temperature()
361 + fluxVarsCache.heatConductionTijFacet()*facetVolVars.temperature();
366 template<
class Problem,
class ElementVolumeVariables >
367 static typename Cache::HeatConductionTransmissibilityContainer
369 const Element& element,
370 const FVElementGeometry& fvGeometry,
371 const ElementVolumeVariables& elemVolVars,
372 const SubControlVolumeFace& scvf)
374 typename Cache::HeatConductionTransmissibilityContainer tij;
375 if (!problem.couplingManager().isCoupled(element, scvf))
378 tij[Cache::insideTijIdx] = ParentType::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf);
383 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
388 if (Dune::FloatCmp::ne(xi, 1.0, 1e-6))
389 DUNE_THROW(Dune::InvalidStateException,
"Xi != 1.0 cannot be used on surface grids");
391 const auto insideScvIdx = scvf.insideScvIdx();
392 const auto& insideScv = fvGeometry.scv(insideScvIdx);
393 const auto& insideVolVars = elemVolVars[insideScvIdx];
394 const auto wIn = Extrusion::area(scvf)
396 insideVolVars.effectiveThermalConductivity(),
397 insideVolVars.extrusionFactor());
400 const auto iBcTypes = problem.interiorBoundaryTypes(element, scvf);
403 if (iBcTypes.hasOnlyNeumann())
408 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
409 const auto wFacet = 2.0*Extrusion::area(scvf)*insideVolVars.extrusionFactor()
410 /sqrt(facetVolVars.extrusionFactor())
411 *
vtmv(scvf.unitOuterNormal(),
412 facetVolVars.effectiveThermalConductivity(),
413 scvf.unitOuterNormal());
415 tij[Cache::insideTijIdx] = wFacet*wIn/(wIn+wFacet);
416 tij[Cache::facetTijIdx] = -tij[Cache::insideTijIdx];
418 else if (iBcTypes.hasOnlyDirichlet())
420 tij[Cache::insideTijIdx] = wIn;
421 tij[Cache::facetTijIdx] = -wIn;
424 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.
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
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:849
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
TpfaFouriersLawCache Cache
export the type for the corresponding cache
Definition: flux/cctpfa/fourierslaw.hh:109
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:177
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:144
Specialization of CCTpfaFacetCouplingFouriersLawImpl for dim<dimWorld.
Definition: multidomain/facet/cellcentered/tpfa/fourierslaw.hh:274
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:345
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:368
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.