24#ifndef DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FOURIERS_LAW_HH
30#include <dune/common/exceptions.hh>
31#include <dune/common/fvector.hh>
32#include <dune/common/float_cmp.hh>
50template<
class TypeTag>
58 using FVElementGeometry =
typename GridGeometry::LocalView;
59 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
60 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
61 using GridView =
typename GridGeometry::GridView;
62 using Element =
typename GridView::template Codim<0>::Entity;
63 using CoordScalar =
typename GridView::ctype;
64 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
66 static constexpr int dim = GridView::dimension;
67 static constexpr int dimWorld = GridView::dimensionworld;
74 template<
class Problem,
class ElementVolumeVariables,
class ElementFluxVarsCache>
75 static Scalar
flux(
const Problem& problem,
76 const Element& element,
77 const FVElementGeometry& fvGeometry,
78 const ElementVolumeVariables& elemVolVars,
79 const SubControlVolumeFace& scvf,
80 const ElementFluxVarsCache& elemFluxVarCache)
83 if (!scvf.interiorBoundary())
84 return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarCache);
86 static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(),
"FacetCoupling.Xi", 1.0);
87 if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) )
88 DUNE_THROW(Dune::NotImplemented,
"Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme");
91 const auto& fluxVarCache = elemFluxVarCache[scvf];
92 const auto& shapeValues = fluxVarCache.shapeValues();
93 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
94 const auto& insideVolVars = elemVolVars[insideScv];
97 const auto bcTypes = problem.interiorBoundaryTypes(element, scvf);
98 if (bcTypes.hasOnlyNeumann())
101 const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf);
105 for (
const auto& scv : scvs(fvGeometry))
106 T += elemVolVars[scv].temperature()*shapeValues[scv.indexInElement()][0];
112 const auto a = facetVolVars.extrusionFactor();
113 auto gradT = scvf.unitOuterNormal();
114 gradT *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a);
115 gradT /= gradT.two_norm2();
116 gradT *= (facetVolVars.temperature() - T);
118 return -1.0*scvf.area()
119 *insideVolVars.extrusionFactor()
120 *
vtmv(scvf.unitOuterNormal(),
121 facetVolVars.effectiveThermalConductivity(),
126 else if (bcTypes.hasOnlyDirichlet())
129 std::vector<Scalar> temperatures(element.subEntities(dim));
130 for (
const auto& scv : scvs(fvGeometry))
131 temperatures[scv.localDofIndex()] = elemVolVars[scv].temperature();
134 for (
const auto& scvfJ : scvfs(fvGeometry))
135 if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement())
136 temperatures[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ]
137 = problem.couplingManager().getLowDimVolVars(element, scvfJ).temperature();
140 Dune::FieldVector<Scalar, dimWorld> gradT(0.0);
141 for (
const auto& scv : scvs(fvGeometry))
142 gradT.axpy(temperatures[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement()));
145 return -1.0*scvf.area()
146 *insideVolVars.extrusionFactor()
147 *
vtmv(scvf.unitOuterNormal(),
148 insideVolVars.effectiveThermalConductivity(),
154 DUNE_THROW(Dune::NotImplemented,
"Mixed boundary types are not supported");
158 template<
class Problem,
class ElementVolumeVariables,
class FluxVarCache>
160 const Element& element,
161 const FVElementGeometry& fvGeometry,
162 const ElementVolumeVariables& elemVolVars,
163 const SubControlVolumeFace& scvf,
164 const FluxVarCache& fluxVarCache,
165 unsigned int phaseIdx)
167 DUNE_THROW(Dune::NotImplemented,
"transmissibilty computation for BoxFacetCouplingFouriersLaw");
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
Specialization of Fourier's Law for the box method.
Definition: flux/box/fourierslaw.hh:45
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: flux/box/fourierslaw.hh:56
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:53
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:75
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:159
Declares all properties used in Dumux.
This file contains the data which is required to calculate energy fluxes due to molecular diffusion w...