version 3.10-dev
box/fourierslawnonequilibrium.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_NONEQUILIBRIUM_HH
14#define DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_NONEQUILIBRIUM_HH
15
16#include <dune/common/fvector.hh>
17
18#include <dumux/common/math.hh>
20
23
25
26namespace Dumux {
27
28// forward declaration
29template <class TypeTag, class DiscretizationMethod>
31
36template <class TypeTag>
37class FouriersLawNonEquilibriumImplementation<TypeTag, DiscretizationMethods::Box>
38{
42 using FVElementGeometry = typename GridGeometry::LocalView;
43 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
44 using Extrusion = Extrusion_t<GridGeometry>;
45 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
46 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
50 using Element = typename GridView::template Codim<0>::Entity;
51
52 static constexpr auto numEnergyEqSolid = getPropValue<TypeTag, Properties::NumEnergyEqSolid>();
53 static constexpr auto numEnergyEqFluid = getPropValue<TypeTag, Properties::NumEnergyEqFluid>();
54 static constexpr auto numEnergyEq = numEnergyEqSolid + numEnergyEqFluid;
55 static constexpr auto sPhaseIdx = ModelTraits::numFluidPhases();
56
57public:
62 static Scalar flux(const Problem& problem,
63 const Element& element,
64 const FVElementGeometry& fvGeometry,
65 const ElementVolumeVariables& elemVolVars,
66 const SubControlVolumeFace& scvf,
67 const int phaseIdx,
68 const ElementFluxVariablesCache& elemFluxVarsCache)
69 {
70 // get inside and outside diffusion tensors and calculate the harmonic mean
71 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
72 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
73 const auto& insideVolVars = elemVolVars[insideScv];
74 const auto& outsideVolVars = elemVolVars[outsideScv];
75 const auto computeLambda = [&](const auto& v){
76 if constexpr (numEnergyEq == 1)
77 return v.effectiveThermalConductivity();
78 else if constexpr (numEnergyEqFluid == 1)
79 return (phaseIdx != sPhaseIdx)
80 ? v.effectiveFluidThermalConductivity()
81 : v.effectiveSolidThermalConductivity();
82 else
83 return v.effectivePhaseThermalConductivity(phaseIdx);
84 };
85
86 auto insideLambda = computeLambda(insideVolVars);
87 auto outsideLambda = computeLambda(outsideVolVars);
88
89 // scale by extrusion factor
90 insideLambda *= insideVolVars.extrusionFactor();
91 outsideLambda *= outsideVolVars.extrusionFactor();
92
93 // the resulting averaged diffusion tensor
94 const auto lambda = faceTensorAverage(insideLambda, outsideLambda, scvf.unitOuterNormal());
95
96 // evaluate gradTemp at integration point
97 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
98
99 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
100 for (auto&& scv : scvs(fvGeometry))
101 {
102 // compute the temperature gradient with the shape functions
103 if (phaseIdx < numEnergyEqFluid)
104 gradTemp.axpy(elemVolVars[scv].temperatureFluid(phaseIdx), fluxVarsCache.gradN(scv.indexInElement()));
105 else
106 gradTemp.axpy(elemVolVars[scv].temperatureSolid(), fluxVarsCache.gradN(scv.indexInElement()));
107 }
108
109 // compute the heat conduction flux
110 return -1.0*vtmv(scvf.unitOuterNormal(), lambda, gradTemp)*Extrusion::area(fvGeometry, scvf);
111 }
112};
113
114} // end namespace Dumux
115
116#endif
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the heat flux within a fluid or solid phase (in J/s) across the given sub-control volume face...
Definition: box/fourierslawnonequilibrium.hh:62
Definition: box/fourierslawnonequilibrium.hh:30
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:880
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
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