3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
ccmpfa/fickslaw.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_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_FICKS_LAW_HH
26
27#include <dumux/common/math.hh>
30
32
33namespace Dumux {
34
36template<class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
37class FicksLawImplementation;
38
43template <class TypeTag, ReferenceSystemFormulation referenceSystem>
44class FicksLawImplementation<TypeTag, DiscretizationMethod::ccmpfa, referenceSystem>
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;
57 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
58 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
59 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
62
63 static constexpr int numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
64 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
65
67 class MpfaFicksLawCacheFiller
68 {
69 public:
72 template<class FluxVariablesCacheFiller>
73 static void fill(FluxVariablesCache& scvfFluxVarsCache,
74 unsigned int phaseIdx, unsigned int compIdx,
75 const Problem& problem,
76 const Element& element,
77 const FVElementGeometry& fvGeometry,
78 const ElementVolumeVariables& elemVolVars,
79 const SubControlVolumeFace& scvf,
80 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
81 {
82 // get interaction volume related data from the filler class & upate the cache
83 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
84 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.secondaryInteractionVolume(),
85 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
86 fluxVarsCacheFiller.secondaryIvDataHandle(),
87 phaseIdx, compIdx);
88 else
89 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.primaryInteractionVolume(),
90 fluxVarsCacheFiller.primaryIvLocalFaceData(),
91 fluxVarsCacheFiller.primaryIvDataHandle(),
92 phaseIdx, compIdx);
93 }
94 };
95
97 class MpfaFicksLawCache
98 {
100 using Stencil = typename DualGridNodalIndexSet::NodalGridStencilType;
101
102 static constexpr int numPhases = GetPropType<TypeTag, Properties::ModelTraits>::numFluidPhases();
103 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
104 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle::DiffusionHandle;
105 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle::DiffusionHandle;
106
108 template< bool doSecondary = considerSecondaryIVs, std::enable_if_t<doSecondary, int> = 0 >
109 void setHandlePointer_(const SecondaryDataHandle& dataHandle)
110 { secondaryHandlePtr_ = &dataHandle; }
111
113 void setHandlePointer_(const PrimaryDataHandle& dataHandle)
114 { primaryHandlePtr_ = &dataHandle; }
115
116 public:
117 // export filler type
118 using Filler = MpfaFicksLawCacheFiller;
119
128 template<class IV, class LocalFaceData, class DataHandle>
129 void updateDiffusion(const IV& iv,
130 const LocalFaceData& localFaceData,
131 const DataHandle& dataHandle,
132 unsigned int phaseIdx, unsigned int compIdx)
133 {
134 stencil_[phaseIdx][compIdx] = &iv.stencil();
135 switchFluxSign_[phaseIdx][compIdx] = localFaceData.isOutsideFace();
136 setHandlePointer_(dataHandle.diffusionHandle());
137 }
138
140 const Stencil& diffusionStencil(unsigned int phaseIdx, unsigned int compIdx) const
141 { return *stencil_[phaseIdx][compIdx]; }
142
144 const PrimaryDataHandle& diffusionPrimaryDataHandle() const { return *primaryHandlePtr_; }
145 const SecondaryDataHandle& diffusionSecondaryDataHandle() const { return *secondaryHandlePtr_; }
146
148 bool diffusionSwitchFluxSign(unsigned int phaseIdx, unsigned int compIdx) const
149 { return switchFluxSign_[phaseIdx][compIdx]; }
150
151
152 private:
154 std::array< std::array<bool, numComponents>, numPhases > switchFluxSign_;
155 std::array< std::array<const Stencil*, numComponents>, numPhases > stencil_;
156
158 const PrimaryDataHandle* primaryHandlePtr_;
159 const SecondaryDataHandle* secondaryHandlePtr_;
160 };
161
162public:
163 // state the discretization method this implementation belongs to
165 //return the reference system
167 { return referenceSystem; }
168 // state the type for the corresponding cache and its filler
169 using Cache = MpfaFicksLawCache;
170
172 static ComponentFluxVector flux(const Problem& problem,
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const SubControlVolumeFace& scvf,
177 const int phaseIdx,
178 const ElementFluxVariablesCache& elemFluxVarsCache)
179 {
180 // obtain this scvf's cache
181 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
182
183 ComponentFluxVector componentFlux(0.0);
184 for (int compIdx = 0; compIdx < numComponents; compIdx++)
185 {
186 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
187 continue;
188
189 // get the scaling factor for the effective diffusive fluxes
190 const auto effFactor = computeEffectivityFactor(elemVolVars, scvf, phaseIdx);
191
192 // if factor is zero, the flux will end up zero anyway
193 if (effFactor == 0.0)
194 continue;
195
196 // calculate the density at the interface
197 const auto rho = interpolateDensity(elemVolVars, scvf, phaseIdx);
198
199 // compute the flux
200 if (fluxVarsCache.usesSecondaryIv())
201 componentFlux[compIdx] = rho*effFactor*computeVolumeFlux(problem,
202 fluxVarsCache,
203 fluxVarsCache.diffusionSecondaryDataHandle(),
204 phaseIdx, compIdx);
205 else
206 componentFlux[compIdx] = rho*effFactor*computeVolumeFlux(problem,
207 fluxVarsCache,
208 fluxVarsCache.diffusionPrimaryDataHandle(),
209 phaseIdx, compIdx);
210 }
211
212 // accumulate the phase component flux
213 for(int compIdx = 0; compIdx < numComponents; compIdx++)
214 if(compIdx != FluidSystem::getMainComponent(phaseIdx) && BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
215 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
216
217 return componentFlux;
218 }
219
220private:
221 template< class Problem, class FluxVarsCache, class DataHandle >
222 static Scalar computeVolumeFlux(const Problem& problem,
223 const FluxVarsCache& cache,
224 const DataHandle& dataHandle,
225 int phaseIdx, int compIdx)
226 {
227 dataHandle.setPhaseIndex(phaseIdx);
228 dataHandle.setComponentIndex(compIdx);
229
230 const bool switchSign = cache.diffusionSwitchFluxSign(phaseIdx, compIdx);
231
232 const auto localFaceIdx = cache.ivLocalFaceIndex();
233 const auto idxInOutside = cache.indexInOutsideFaces();
234 const auto& xj = dataHandle.uj();
235 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
236 : (!switchSign ? dataHandle.T()[localFaceIdx]
237 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
238 Scalar scvfFlux = tij*xj;
239
240 // switch the sign if necessary
241 if (switchSign)
242 scvfFlux *= -1.0;
243
244 return scvfFlux;
245 }
246
248 static Scalar interpolateDensity(const ElementVolumeVariables& elemVolVars,
249 const SubControlVolumeFace& scvf,
250 const unsigned int phaseIdx)
251 {
252 // use arithmetic mean of the densities around the scvf
253 if (!scvf.boundary())
254 {
255 const Scalar rhoInside = massOrMolarDensity(elemVolVars[scvf.insideScvIdx()], referenceSystem, phaseIdx);
256
257 Scalar rho = rhoInside;
258 for (const auto outsideIdx : scvf.outsideScvIndices())
259 {
260 const Scalar rhoOutside = massOrMolarDensity(elemVolVars[outsideIdx], referenceSystem, phaseIdx);
261 rho += rhoOutside;
262 }
263 return rho/(scvf.outsideScvIndices().size()+1);
264 }
265 else
266 return massOrMolarDensity(elemVolVars[scvf.outsideScvIdx()], referenceSystem, phaseIdx);
267 }
268
273 static Scalar computeEffectivityFactor(const ElementVolumeVariables& elemVolVars,
274 const SubControlVolumeFace& scvf,
275 const unsigned int phaseIdx)
276 {
277 using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
278
279 // use the harmonic mean between inside and outside
280 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
281 const auto factor = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
282 insideVolVars.saturation(phaseIdx),
283 /*Diffusion coefficient*/ 1.0);
284
285 if (!scvf.boundary())
286 {
287 // interpret outside factor as arithmetic mean
288 Scalar outsideFactor = 0.0;
289 for (const auto outsideIdx : scvf.outsideScvIndices())
290 {
291 const auto& outsideVolVars = elemVolVars[outsideIdx];
292 outsideFactor += EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
293 outsideVolVars.saturation(phaseIdx),
294 /*Diffusion coefficient*/ 1.0);
295 }
296 outsideFactor /= scvf.outsideScvIndices().size();
297
298 // use the harmonic mean of the two
299 return harmonicMean(factor, outsideFactor);
300 }
301
302 return factor;
303 }
304};
305
306} // end namespace
307
308#endif
Define some often used mathematical functions.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:55
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
make the local view function available whenever we use the grid geometry
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
forward declaration of the method-specific implemetation
Definition: box/fickslaw.hh:38
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: ccmpfa/fickslaw.hh:166
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Compute the diffusive flux across an scvf.
Definition: ccmpfa/fickslaw.hh:172
Declares all properties used in Dumux.