3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cctpfa/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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_NONEQUILIBRIUM_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FOURIERS_LAW_NONEQUILIBRIUM_HH
26
32
33namespace Dumux {
34
35// forward declaration
36template<class TypeTag, DiscretizationMethod discMethod>
37class FouriersLawNonEquilibriumImplementation;
38
43template <class TypeTag>
45{
50 using FVElementGeometry = typename GridGeometry::LocalView;
51 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
52 using Extrusion = Extrusion_t<GridGeometry>;
53 using GridView = typename GridGeometry::GridView;
54 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
55 using Element = typename GridView::template Codim<0>::Entity;
57
58 static constexpr int dim = GridView::dimension;
59 static constexpr int dimWorld = GridView::dimensionworld;
60
62
63 static constexpr auto numEnergyEqSolid = getPropValue<TypeTag, Properties::NumEnergyEqSolid>();
64 static constexpr auto numEnergyEqFluid = getPropValue<TypeTag, Properties::NumEnergyEqFluid>();
65 static constexpr auto numEnergyEq = numEnergyEqSolid + numEnergyEqFluid;
66 static constexpr auto sPhaseIdx = ModelTraits::numFluidPhases();
67
68public:
71
73
78 static Scalar flux(const Problem& problem,
79 const Element& element,
80 const FVElementGeometry& fvGeometry,
81 const ElementVolumeVariables& elemVolVars,
82 const SubControlVolumeFace& scvf,
83 const int phaseIdx,
84 const ElementFluxVarsCache& elemFluxVarsCache)
85 {
86 Scalar tInside = 0.0;
87 Scalar tOutside = 0.0;
88 // get the inside/outside temperatures
89 if (phaseIdx < numEnergyEqFluid)
90 {
91 tInside += elemVolVars[scvf.insideScvIdx()].temperatureFluid(phaseIdx);
92 tOutside += elemVolVars[scvf.outsideScvIdx()].temperatureFluid(phaseIdx);
93 }
94 else //temp solid
95 {
96 tInside += elemVolVars[scvf.insideScvIdx()].temperatureSolid();
97 tOutside += elemVolVars[scvf.outsideScvIdx()].temperatureSolid();
98 }
99
100 Scalar tij = calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx);
101 return tij*(tInside - tOutside);
102 }
103
105 static Scalar calculateTransmissibility(const Problem& problem,
106 const Element& element,
107 const FVElementGeometry& fvGeometry,
108 const ElementVolumeVariables& elemVolVars,
109 const SubControlVolumeFace& scvf,
110 const int phaseIdx)
111 {
112 const auto insideScvIdx = scvf.insideScvIdx();
113 const auto& insideScv = fvGeometry.scv(insideScvIdx);
114 const auto& insideVolVars = elemVolVars[insideScvIdx];
115 const auto computeLambda = [&](const auto& v){
116 if constexpr (numEnergyEq == 1)
117 return v.effectiveThermalConductivity();
118 else if constexpr (numEnergyEqFluid == 1)
119 return (phaseIdx != sPhaseIdx)
120 ? v.effectiveFluidThermalConductivity()
121 : v.effectiveSolidThermalConductivity();
122 else
123 return v.effectivePhaseThermalConductivity(phaseIdx);
124 };
125
126 const auto insideLambda = computeLambda(insideVolVars);
127 const Scalar ti = computeTpfaTransmissibility(scvf, insideScv, insideLambda, insideVolVars.extrusionFactor());
128
129 // for the boundary (dirichlet) or at branching points we only need ti
130 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
131 return Extrusion::area(scvf)*ti;
132 else // otherwise we compute a tpfa harmonic mean
133 {
134 const auto outsideScvIdx = scvf.outsideScvIdx();
135 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
136 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
137 const auto outsideLambda = computeLambda(outsideVolVars);
138
139 Scalar tj;
140 if (dim == dimWorld)
141 // assume the normal vector from outside is anti parallel so we save flipping a vector
142 tj = -1.0*computeTpfaTransmissibility(scvf, outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
143 else
144 tj = computeTpfaTransmissibility(fvGeometry.flipScvf(scvf.index()), outsideScv, outsideLambda, outsideVolVars.extrusionFactor());
145
146 // check for division by zero!
147 if (ti*tj <= 0.0)
148 return 0.0;
149 else
150 return Extrusion::area(scvf)*(ti * tj)/(ti + tj);
151 }
152 }
153};
154
155} // end namespace Dumux
156
157#endif
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
Classes related to flux variables caching.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:47
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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
Fourier's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: cctpfa/fourierslawnonequilibrium.hh:45
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx)
Compute transmissibilities.
Definition: cctpfa/fourierslawnonequilibrium.hh:105
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVarsCache &elemFluxVarsCache)
Returns the heat flux within a fluid or solid phase (in J/s) across the given sub-control volume face...
Definition: cctpfa/fourierslawnonequilibrium.hh:78
Definition: fluxvariablescaching.hh:65
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...