version 3.8
flux/box/fourierslaw.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_HH
14#define DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_HH
15
16#include <dumux/common/math.hh>
21
22namespace Dumux {
23
24// forward declaration
25template<class TypeTag, class DiscretizationMethod>
27
32template <class TypeTag>
33class FouriersLawImplementation<TypeTag, DiscretizationMethods::Box>
34{
38 using FVElementGeometry = typename GridGeometry::LocalView;
39 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
40 using Extrusion = Extrusion_t<GridGeometry>;
41 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
42 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
44 using Element = typename GridView::template Codim<0>::Entity;
45
46public:
54 static Scalar flux(const Problem& problem,
55 const Element& element,
56 const FVElementGeometry& fvGeometry,
57 const ElementVolumeVariables& elemVolVars,
58 const SubControlVolumeFace& scvf,
59 const ElementFluxVariablesCache& elemFluxVarsCache)
60 {
61 // get inside and outside diffusion tensors and calculate the harmonic mean
62 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
63 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
64 const auto& insideVolVars = elemVolVars[insideScv];
65 const auto& outsideVolVars = elemVolVars[outsideScv];
66
67 // effective diffusion tensors
68 auto insideLambda = insideVolVars.effectiveThermalConductivity();
69 auto outsideLambda = outsideVolVars.effectiveThermalConductivity();
70
71 // scale by extrusion factor
72 insideLambda *= insideVolVars.extrusionFactor();
73 outsideLambda *= outsideVolVars.extrusionFactor();
74
75 // the resulting averaged diffusion tensor
76 const auto lambda = faceTensorAverage(insideLambda, outsideLambda, scvf.unitOuterNormal());
77
78 // evaluate gradTemp at integration point
79 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
80
81 // compute the temperature gradient with the shape functions
82 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
83 for (auto&& scv : scvs(fvGeometry))
84 gradTemp.axpy(elemVolVars[scv].temperature(), fluxVarsCache.gradN(scv.indexInElement()));
85
86 // compute the heat conduction flux
87 return -1.0*vtmv(scvf.unitOuterNormal(), lambda, gradTemp)*Extrusion::area(fvGeometry, scvf);
88 }
89};
90
91} // end namespace Dumux
92
93#endif
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.
A free function to average a Tensor at an interface.
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.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:39
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:29