3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_NONEQUILIBRIUM_HH
26#define DUMUX_DISCRETIZATION_BOX_FOURIERS_LAW_NONEQUILIBRIUM_HH
27
28#include <dune/common/fvector.hh>
29
30#include <dumux/common/math.hh>
33
35
36
37namespace Dumux {
38
39// forward declaration
40template <class TypeTag, DiscretizationMethod discMethod>
41class FouriersLawNonEquilibriumImplementation;
42
47template <class TypeTag>
49{
52 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
53 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
55 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
59 using Element = typename GridView::template Codim<0>::Entity;
60
61 static constexpr auto numEnergyEqSolid = getPropValue<TypeTag, Properties::NumEnergyEqSolid>();
62 static constexpr auto numEnergyEqFluid = getPropValue<TypeTag, Properties::NumEnergyEqFluid>();
63 static constexpr auto numEnergyEq = numEnergyEqSolid + numEnergyEqFluid;
64 static constexpr auto sPhaseIdx = ModelTraits::numFluidPhases();
65
66public:
67 static Scalar flux(const Problem& problem,
68 const Element& element,
69 const FVElementGeometry& fvGeometry,
70 const ElementVolumeVariables& elemVolVars,
71 const SubControlVolumeFace& scvf,
72 const int phaseIdx,
73 const ElementFluxVariablesCache& elemFluxVarsCache)
74 {
75 // get inside and outside diffusion tensors and calculate the harmonic mean
76 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
77 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
78 const auto& insideVolVars = elemVolVars[insideScv];
79 const auto& outsideVolVars = elemVolVars[outsideScv];
80 const auto computeLambda = [&](const auto& v){
81 if constexpr (numEnergyEq == 1)
82 return v.effectiveThermalConductivity();
83 else if constexpr (numEnergyEqFluid == 1)
84 return (phaseIdx != sPhaseIdx)
85 ? v.effectiveFluidThermalConductivity()
86 : v.effectiveSolidThermalConductivity();
87 else
88 return v.effectivePhaseThermalConductivity(phaseIdx);
89 };
90
91 auto insideLambda = computeLambda(insideVolVars);
92 auto outsideLambda = computeLambda(outsideVolVars);
93
94 // scale by extrusion factor
95 insideLambda *= insideVolVars.extrusionFactor();
96 outsideLambda *= outsideVolVars.extrusionFactor();
97
98 // the resulting averaged diffusion tensor
99 const auto lambda = problem.spatialParams().harmonicMean(insideLambda, outsideLambda, scvf.unitOuterNormal());
100
101 // evaluate gradTemp at integration point
102 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
103
104 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
105 for (auto&& scv : scvs(fvGeometry))
106 {
107 // compute the temperature gradient with the shape functions
108 if (phaseIdx < numEnergyEqFluid)
109 gradTemp.axpy(elemVolVars[scv].temperatureFluid(phaseIdx), fluxVarsCache.gradN(scv.indexInElement()));
110 else
111 gradTemp.axpy(elemVolVars[scv].temperatureSolid(), fluxVarsCache.gradN(scv.indexInElement()));
112 }
113
114 // comute the heat conduction flux
115 return -1.0*vtmv(scvf.unitOuterNormal(), lambda, gradTemp)*scvf.area();
116 }
117};
118
119} // end namespace Dumux
120
121#endif
Define some often used mathematical functions.
Helpers for deprecation.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
Definition: adapt.hh:29
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
Definition: fourierslawnonequilibrium.hh:36
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: box/fourierslawnonequilibrium.hh:67
Declares all properties used in Dumux.