3.6-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>
32
35
37
38namespace Dumux {
39
40// forward declaration
41template <class TypeTag, class DiscretizationMethod>
43
48template <class TypeTag>
49class FouriersLawNonEquilibriumImplementation<TypeTag, DiscretizationMethods::Box>
50{
54 using FVElementGeometry = typename GridGeometry::LocalView;
55 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
56 using Extrusion = Extrusion_t<GridGeometry>;
57 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
58 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
62 using Element = typename GridView::template Codim<0>::Entity;
63
64 static constexpr auto numEnergyEqSolid = getPropValue<TypeTag, Properties::NumEnergyEqSolid>();
65 static constexpr auto numEnergyEqFluid = getPropValue<TypeTag, Properties::NumEnergyEqFluid>();
66 static constexpr auto numEnergyEq = numEnergyEqSolid + numEnergyEqFluid;
67 static constexpr auto sPhaseIdx = ModelTraits::numFluidPhases();
68
69public:
74 static Scalar flux(const Problem& problem,
75 const Element& element,
76 const FVElementGeometry& fvGeometry,
77 const ElementVolumeVariables& elemVolVars,
78 const SubControlVolumeFace& scvf,
79 const int phaseIdx,
80 const ElementFluxVariablesCache& elemFluxVarsCache)
81 {
82 // get inside and outside diffusion tensors and calculate the harmonic mean
83 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
84 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
85 const auto& insideVolVars = elemVolVars[insideScv];
86 const auto& outsideVolVars = elemVolVars[outsideScv];
87 const auto computeLambda = [&](const auto& v){
88 if constexpr (numEnergyEq == 1)
89 return v.effectiveThermalConductivity();
90 else if constexpr (numEnergyEqFluid == 1)
91 return (phaseIdx != sPhaseIdx)
92 ? v.effectiveFluidThermalConductivity()
93 : v.effectiveSolidThermalConductivity();
94 else
95 return v.effectivePhaseThermalConductivity(phaseIdx);
96 };
97
98 auto insideLambda = computeLambda(insideVolVars);
99 auto outsideLambda = computeLambda(outsideVolVars);
100
101 // scale by extrusion factor
102 insideLambda *= insideVolVars.extrusionFactor();
103 outsideLambda *= outsideVolVars.extrusionFactor();
104
105 // the resulting averaged diffusion tensor
106 const auto lambda = faceTensorAverage(insideLambda, outsideLambda, scvf.unitOuterNormal());
107
108 // evaluate gradTemp at integration point
109 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
110
111 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
112 for (auto&& scv : scvs(fvGeometry))
113 {
114 // compute the temperature gradient with the shape functions
115 if (phaseIdx < numEnergyEqFluid)
116 gradTemp.axpy(elemVolVars[scv].temperatureFluid(phaseIdx), fluxVarsCache.gradN(scv.indexInElement()));
117 else
118 gradTemp.axpy(elemVolVars[scv].temperatureSolid(), fluxVarsCache.gradN(scv.indexInElement()));
119 }
120
121 // compute the heat conduction flux
122 return -1.0*vtmv(scvf.unitOuterNormal(), lambda, gradTemp)*Extrusion::area(fvGeometry, scvf);
123 }
124};
125
126} // end namespace Dumux
127
128#endif
A free function to average a Tensor at an interface.
Define some often used mathematical functions.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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:863
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
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:41
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
Definition: box/fourierslawnonequilibrium.hh:42
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:74
Declares all properties used in Dumux.