3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
flux/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 <dune/common/fvector.hh>
30
33
34namespace Dumux {
35
37template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
38class FicksLawImplementation;
39
44template <class TypeTag, ReferenceSystemFormulation referenceSystem>
45class FicksLawImplementation<TypeTag, DiscretizationMethods::CCMpfa, referenceSystem>
46{
50 using Element = typename GridView::template Codim<0>::Entity;
51
52 static constexpr int dim = GridView::dimension;
53 static constexpr int dimWorld = GridView::dimensionworld;
54
56 using FVElementGeometry = typename GridGeometry::LocalView;
58 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
59 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
61 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
62 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
64
65 static constexpr int numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
67 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
68
70 class MpfaFicksLawCacheFiller
71 {
72 public:
75 template<class FluxVariablesCacheFiller>
76 static void fill(FluxVariablesCache& scvfFluxVarsCache,
77 unsigned int phaseIdx, unsigned int compIdx,
78 const Problem& problem,
79 const Element& element,
80 const FVElementGeometry& fvGeometry,
81 const ElementVolumeVariables& elemVolVars,
82 const SubControlVolumeFace& scvf,
83 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
84 {
85 // get interaction volume related data from the filler class & update the cache
86 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
87 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.secondaryInteractionVolume(),
88 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
89 fluxVarsCacheFiller.secondaryIvDataHandle(),
90 phaseIdx, compIdx);
91 else
92 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.primaryInteractionVolume(),
93 fluxVarsCacheFiller.primaryIvLocalFaceData(),
94 fluxVarsCacheFiller.primaryIvDataHandle(),
95 phaseIdx, compIdx);
96 }
97 };
98
100 class MpfaFicksLawCache
101 {
103 using Stencil = typename DualGridNodalIndexSet::NodalGridStencilType;
104
105 static constexpr int numPhases = GetPropType<TypeTag, Properties::ModelTraits>::numFluidPhases();
106 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
107 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle::DiffusionHandle;
108 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle::DiffusionHandle;
109
111 template< bool doSecondary = considerSecondaryIVs, std::enable_if_t<doSecondary, int> = 0 >
112 void setHandlePointer_(const SecondaryDataHandle& dataHandle)
113 { secondaryHandlePtr_ = &dataHandle; }
114
116 void setHandlePointer_(const PrimaryDataHandle& dataHandle)
117 { primaryHandlePtr_ = &dataHandle; }
118
119 public:
120 // export filler type
121 using Filler = MpfaFicksLawCacheFiller;
122
131 template<class IV, class LocalFaceData, class DataHandle>
132 void updateDiffusion(const IV& iv,
133 const LocalFaceData& localFaceData,
134 const DataHandle& dataHandle,
135 unsigned int phaseIdx, unsigned int compIdx)
136 {
137 stencil_[phaseIdx][compIdx] = &iv.stencil();
138 switchFluxSign_[phaseIdx][compIdx] = localFaceData.isOutsideFace();
139 setHandlePointer_(dataHandle.diffusionHandle());
140 }
141
143 const Stencil& diffusionStencil(unsigned int phaseIdx, unsigned int compIdx) const
144 { return *stencil_[phaseIdx][compIdx]; }
145
147 const PrimaryDataHandle& diffusionPrimaryDataHandle() const { return *primaryHandlePtr_; }
148 const SecondaryDataHandle& diffusionSecondaryDataHandle() const { return *secondaryHandlePtr_; }
149
151 bool diffusionSwitchFluxSign(unsigned int phaseIdx, unsigned int compIdx) const
152 { return switchFluxSign_[phaseIdx][compIdx]; }
153
154
155 private:
157 std::array< std::array<bool, numComponents>, numPhases > switchFluxSign_;
158 std::array< std::array<const Stencil*, numComponents>, numPhases > stencil_;
159
161 const PrimaryDataHandle* primaryHandlePtr_;
162 const SecondaryDataHandle* secondaryHandlePtr_;
163 };
164
165public:
167 // state the discretization method this implementation belongs to
168 static constexpr DiscretizationMethod discMethod{};
169
170 //return the reference system
172 { return referenceSystem; }
173
174 // state the type for the corresponding cache and its filler
175 using Cache = MpfaFicksLawCache;
176
177 // export the diffusion container
179
186 static ComponentFluxVector flux(const Problem& problem,
187 const Element& element,
188 const FVElementGeometry& fvGeometry,
189 const ElementVolumeVariables& elemVolVars,
190 const SubControlVolumeFace& scvf,
191 const int phaseIdx,
192 const ElementFluxVariablesCache& elemFluxVarsCache)
193 {
194 // obtain this scvf's cache
195 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
196
197 ComponentFluxVector componentFlux(0.0);
198 for (int compIdx = 0; compIdx < numComponents; compIdx++)
199 {
200 if constexpr (!FluidSystem::isTracerFluidSystem())
201 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
202 continue;
203
204 // calculate the density at the interface
205 const auto rho = interpolateDensity(elemVolVars, scvf, phaseIdx);
206
207 // compute the flux
208 if (fluxVarsCache.usesSecondaryIv())
209 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
210 fluxVarsCache,
211 fluxVarsCache.diffusionSecondaryDataHandle(),
212 phaseIdx, compIdx);
213 else
214 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
215 fluxVarsCache,
216 fluxVarsCache.diffusionPrimaryDataHandle(),
217 phaseIdx, compIdx);
218 }
219
220 // accumulate the phase component flux
221 for (int compIdx = 0; compIdx < numComponents; compIdx++)
222 if constexpr (!FluidSystem::isTracerFluidSystem())
223 if (compIdx != FluidSystem::getMainComponent(phaseIdx) && BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
224 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
225
226 return componentFlux;
227 }
228
229private:
230 template< class Problem, class FluxVarsCache, class DataHandle >
231 static Scalar computeVolumeFlux(const Problem& problem,
232 const FluxVarsCache& cache,
233 const DataHandle& dataHandle,
234 int phaseIdx, int compIdx)
235 {
236 dataHandle.setPhaseIndex(phaseIdx);
237 dataHandle.setComponentIndex(compIdx);
238
239 const bool switchSign = cache.diffusionSwitchFluxSign(phaseIdx, compIdx);
240
241 const auto localFaceIdx = cache.ivLocalFaceIndex();
242 const auto idxInOutside = cache.indexInOutsideFaces();
243 const auto& xj = dataHandle.uj();
244 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
245 : (!switchSign ? dataHandle.T()[localFaceIdx]
246 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
247 Scalar scvfFlux = tij*xj;
248
249 // switch the sign if necessary
250 if (switchSign)
251 scvfFlux *= -1.0;
252
253 return scvfFlux;
254 }
255
257 static Scalar interpolateDensity(const ElementVolumeVariables& elemVolVars,
258 const SubControlVolumeFace& scvf,
259 const unsigned int phaseIdx)
260 {
261 // use arithmetic mean of the densities around the scvf
262 if (!scvf.boundary())
263 {
264 const Scalar rhoInside = massOrMolarDensity(elemVolVars[scvf.insideScvIdx()], referenceSystem, phaseIdx);
265
266 Scalar rho = rhoInside;
267 for (const auto outsideIdx : scvf.outsideScvIndices())
268 {
269 const Scalar rhoOutside = massOrMolarDensity(elemVolVars[outsideIdx], referenceSystem, phaseIdx);
270 rho += rhoOutside;
271 }
272 return rho/(scvf.outsideScvIndices().size()+1);
273 }
274 else
275 return massOrMolarDensity(elemVolVars[scvf.outsideScvIdx()], referenceSystem, phaseIdx);
276 }
277};
278
279} // end namespace
280
281#endif
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
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
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/fickslaw.hh:44
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: flux/ccmpfa/fickslaw.hh:186
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/ccmpfa/fickslaw.hh:171
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:44
Declares all properties used in Dumux.