3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/ccmpfa/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 * 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_MPFA_FOURIERS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_FOURIERS_LAW_HH
26
27#include <dune/common/dynvector.hh>
28#include <dune/common/dynmatrix.hh>
29
32
33namespace Dumux {
34
36template<class TypeTag, class DiscretizationMethod>
37class FouriersLawImplementation;
38
43template <class TypeTag>
44class FouriersLawImplementation<TypeTag, DiscretizationMethods::CCMpfa>
45{
49 using Element = typename GridView::template Codim<0>::Entity;
50
51 static constexpr int dim = GridView::dimension;
52 static constexpr int dimWorld = GridView::dimensionworld;
53
55 using FVElementGeometry = typename GridGeometry::LocalView;
56 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
57 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
59 using ElementFluxVarsCache = typename GridFluxVariablesCache::LocalView;
60 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
61
63 class MpfaFouriersLawCacheFiller
64 {
65 public:
68 template<class FluxVariablesCacheFiller>
69 static void fill(FluxVariablesCache& scvfFluxVarsCache,
70 const Problem& problem,
71 const Element& element,
72 const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars,
74 const SubControlVolumeFace& scvf,
75 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
76 {
77 // get interaction volume from the flux vars cache filler & update the cache
78 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
79 scvfFluxVarsCache.updateHeatConduction(fluxVarsCacheFiller.secondaryInteractionVolume(),
80 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
81 fluxVarsCacheFiller.secondaryIvDataHandle());
82 else
83 scvfFluxVarsCache.updateHeatConduction(fluxVarsCacheFiller.primaryInteractionVolume(),
84 fluxVarsCacheFiller.primaryIvLocalFaceData(),
85 fluxVarsCacheFiller.primaryIvDataHandle());
86 }
87 };
88
90 class MpfaFouriersLawCache
91 {
93 using Stencil = typename DualGridNodalIndexSet::NodalGridStencilType;
94
95 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
96 using PrimaryDataHandle = typename ElementFluxVarsCache::PrimaryIvDataHandle::HeatConductionHandle;
97 using SecondaryDataHandle = typename ElementFluxVarsCache::SecondaryIvDataHandle::HeatConductionHandle;
98
100 template< bool doSecondary = considerSecondaryIVs, std::enable_if_t<doSecondary, int> = 0 >
101 void setHandlePointer_(const SecondaryDataHandle& dataHandle)
102 { secondaryHandlePtr_ = &dataHandle; }
103
105 void setHandlePointer_(const PrimaryDataHandle& dataHandle)
106 { primaryHandlePtr_ = &dataHandle; }
107
108 public:
109 // export filler type
110 using Filler = MpfaFouriersLawCacheFiller;
111
120 template<class IV, class LocalFaceData, class DataHandle>
121 void updateHeatConduction(const IV& iv,
122 const LocalFaceData& localFaceData,
123 const DataHandle& dataHandle)
124 {
125 switchFluxSign_ = localFaceData.isOutsideFace();
126 stencil_ = &iv.stencil();
127 setHandlePointer_(dataHandle.heatConductionHandle());
128 }
129
131 const Stencil& heatConductionStencil() const { return *stencil_; }
132
134 const PrimaryDataHandle& heatConductionPrimaryDataHandle() const { return *primaryHandlePtr_; }
135 const SecondaryDataHandle& heatConductionSecondaryDataHandle() const { return *secondaryHandlePtr_; }
136
138 bool heatConductionSwitchFluxSign() const { return switchFluxSign_; }
139
140 private:
141 bool switchFluxSign_;
142
144 const PrimaryDataHandle* primaryHandlePtr_;
145 const SecondaryDataHandle* secondaryHandlePtr_;
146
148 const Stencil* stencil_;
149 };
150
151public:
153 // state the discretization method this implementation belongs to
154 static constexpr DiscretizationMethod discMethod{};
155
156 // state the type for the corresponding cache and its filler
157 using Cache = MpfaFouriersLawCache;
158
166 static Scalar flux(const Problem& problem,
167 const Element& element,
168 const FVElementGeometry& fvGeometry,
169 const ElementVolumeVariables& elemVolVars,
170 const SubControlVolumeFace& scvf,
171 const ElementFluxVarsCache& elemFluxVarsCache)
172 {
173 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
174
175 // forward to the private function taking the iv data handle
176 if (fluxVarsCache.usesSecondaryIv())
177 return flux_(problem, fluxVarsCache, fluxVarsCache.heatConductionSecondaryDataHandle());
178 else
179 return flux_(problem, fluxVarsCache, fluxVarsCache.heatConductionPrimaryDataHandle());
180 }
181
182private:
183 template< class Problem, class FluxVarsCache, class DataHandle >
184 static Scalar flux_(const Problem& problem,
185 const FluxVarsCache& cache,
186 const DataHandle& dataHandle)
187 {
188 const bool switchSign = cache.heatConductionSwitchFluxSign();
189
190 const auto localFaceIdx = cache.ivLocalFaceIndex();
191 const auto idxInOutside = cache.indexInOutsideFaces();
192 const auto& Tj = dataHandle.uj();
193 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
194 : (!switchSign ? dataHandle.T()[localFaceIdx]
195 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
196 Scalar scvfFlux = tij*Tj;
197
198 // switch the sign if necessary
199 if (switchSign)
200 scvfFlux *= -1.0;
201
202 return scvfFlux;
203 }
204};
205
206} // end namespace Dumux
207
208#endif
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Definition: method.hh:45
forward declaration of the method-specific implementation
Definition: flux/box/fourierslaw.hh:38
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVarsCache &elemFluxVarsCache)
Returns the heat flux within the porous medium (in J/s) across the given sub-control volume face.
Definition: flux/ccmpfa/fourierslaw.hh:166
MpfaFouriersLawCache Cache
Definition: flux/ccmpfa/fourierslaw.hh:157
Declares all properties used in Dumux.