12#ifndef DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FOURIERS_LAW_HH
13#define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FOURIERS_LAW_HH
18#include <dune/common/exceptions.hh>
19#include <dune/common/fvector.hh>
20#include <dune/common/float_cmp.hh>
39template<
class TypeTag>
47 using FVElementGeometry =
typename GridGeometry::LocalView;
48 using SubControlVolume =
typename GridGeometry::SubControlVolume;
49 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
51 using GridView =
typename GridGeometry::GridView;
52 using Element =
typename GridView::template Codim<0>::Entity;
53 using CoordScalar =
typename GridView::ctype;
54 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
56 static constexpr int dim = GridView::dimension;
57 static constexpr int dimWorld = GridView::dimensionworld;
64 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
65 static Scalar
flux(
const Problem& problem,
66 const Element& element,
67 const FVElementGeometry& fvGeometry,
68 const ElementVolumeVariables& elemVolVars,
69 const SubControlVolumeFace& scvf,
70 const ElementFluxVarsCache& elemFluxVarCache)
73 if (!scvf.interiorBoundary())
74 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarCache);
76 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
77 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
78 DUNE_THROW(Dune::NotImplemented,
"Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme");
81 const auto& fluxVarCache = elemFluxVarCache[scvf];
82 const auto& shapeValues = fluxVarCache.shapeValues();
83 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
84 const auto& insideVolVars = elemVolVars[insideScv];
87 const auto bcTypes = problem.interiorBoundaryTypes(element, scvf);
88 if (bcTypes.hasOnlyNeumann())
91 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
95 for (
const auto& scv : scvs(fvGeometry))
96 T += elemVolVars[scv].temperature()*shapeValues[scv.indexInElement()][0];
102 const auto a = facetVolVars.extrusionFactor();
103 auto gradT = scvf.unitOuterNormal();
104 gradT *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a);
105 gradT /= gradT.two_norm2();
106 gradT *= (facetVolVars.temperature() - T);
108 return -1.0*Extrusion::area(fvGeometry, scvf)
109 *insideVolVars.extrusionFactor()
110 *
vtmv(scvf.unitOuterNormal(),
111 facetVolVars.effectiveThermalConductivity(),
116 else if (bcTypes.hasOnlyDirichlet())
119 std::vector<Scalar> temperatures(element.subEntities(dim));
120 for (
const auto& scv : scvs(fvGeometry))
121 temperatures[scv.localDofIndex()] = elemVolVars[scv].temperature();
124 for (
const auto& scvfJ : scvfs(fvGeometry))
125 if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement())
126 temperatures[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ]
127 = problem.couplingManager().getLowDimVolVars(element, scvfJ).temperature();
130 Dune::FieldVector<Scalar, dimWorld> gradT(0.0);
131 for (
const auto& scv : scvs(fvGeometry))
132 gradT.axpy(temperatures[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement()));
135 return -1.0*Extrusion::area(fvGeometry, scvf)
136 *insideVolVars.extrusionFactor()
137 *
vtmv(scvf.unitOuterNormal(),
138 insideVolVars.effectiveThermalConductivity(),
144 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary types are not supported");
148 template<
class Problem,
class ElementVolumeVariables,
class FluxVarCache>
150 const Element& element,
151 const FVElementGeometry& fvGeometry,
152 const ElementVolumeVariables& elemVolVars,
153 const SubControlVolumeFace& scvf,
154 const FluxVarCache& fluxVarCache,
155 unsigned int phaseIdx)
157 DUNE_THROW(Dune::NotImplemented,
"transmissibilty computation for BoxFacetCouplingFouriersLaw");
Fourier's law for the box scheme in the context of coupled models where coupling occurs across the fa...
Definition: multidomain/facet/box/fourierslaw.hh:42
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarCache)
Compute the conductive heat flux across a sub-control volume face.
Definition: multidomain/facet/box/fourierslaw.hh:65
static std::vector< Scalar > calculateTransmissibilities(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const FluxVarCache &fluxVarCache, unsigned int phaseIdx)
Definition: multidomain/facet/box/fourierslaw.hh:149
Specialization of Fourier's Law for the box method.
Definition: flux/box/fourierslaw.hh:34
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the heat flux within the porous medium (in J/s) across the given sub-control volume face.
Definition: flux/box/fourierslaw.hh:54
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.
This file contains the data which is required to calculate energy fluxes due to molecular diffusion w...
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:851
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.