3.2-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, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
38class FicksLawImplementation;
39
44template <class TypeTag, ReferenceSystemFormulation referenceSystem>
45class FicksLawImplementation<TypeTag, DiscretizationMethod::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;
60 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
63
64 static constexpr int numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
66 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
67
69 class MpfaFicksLawCacheFiller
70 {
71 public:
74 template<class FluxVariablesCacheFiller>
75 static void fill(FluxVariablesCache& scvfFluxVarsCache,
76 unsigned int phaseIdx, unsigned int compIdx,
77 const Problem& problem,
78 const Element& element,
79 const FVElementGeometry& fvGeometry,
80 const ElementVolumeVariables& elemVolVars,
81 const SubControlVolumeFace& scvf,
82 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
83 {
84 // get interaction volume related data from the filler class & upate the cache
85 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
86 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.secondaryInteractionVolume(),
87 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
88 fluxVarsCacheFiller.secondaryIvDataHandle(),
89 phaseIdx, compIdx);
90 else
91 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.primaryInteractionVolume(),
92 fluxVarsCacheFiller.primaryIvLocalFaceData(),
93 fluxVarsCacheFiller.primaryIvDataHandle(),
94 phaseIdx, compIdx);
95 }
96 };
97
99 class MpfaFicksLawCache
100 {
102 using Stencil = typename DualGridNodalIndexSet::NodalGridStencilType;
103
104 static constexpr int numPhases = GetPropType<TypeTag, Properties::ModelTraits>::numFluidPhases();
105 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
106 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle::DiffusionHandle;
107 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle::DiffusionHandle;
108
110 template< bool doSecondary = considerSecondaryIVs, std::enable_if_t<doSecondary, int> = 0 >
111 void setHandlePointer_(const SecondaryDataHandle& dataHandle)
112 { secondaryHandlePtr_ = &dataHandle; }
113
115 void setHandlePointer_(const PrimaryDataHandle& dataHandle)
116 { primaryHandlePtr_ = &dataHandle; }
117
118 public:
119 // export filler type
120 using Filler = MpfaFicksLawCacheFiller;
121
130 template<class IV, class LocalFaceData, class DataHandle>
131 void updateDiffusion(const IV& iv,
132 const LocalFaceData& localFaceData,
133 const DataHandle& dataHandle,
134 unsigned int phaseIdx, unsigned int compIdx)
135 {
136 stencil_[phaseIdx][compIdx] = &iv.stencil();
137 switchFluxSign_[phaseIdx][compIdx] = localFaceData.isOutsideFace();
138 setHandlePointer_(dataHandle.diffusionHandle());
139 }
140
142 const Stencil& diffusionStencil(unsigned int phaseIdx, unsigned int compIdx) const
143 { return *stencil_[phaseIdx][compIdx]; }
144
146 const PrimaryDataHandle& diffusionPrimaryDataHandle() const { return *primaryHandlePtr_; }
147 const SecondaryDataHandle& diffusionSecondaryDataHandle() const { return *secondaryHandlePtr_; }
148
150 bool diffusionSwitchFluxSign(unsigned int phaseIdx, unsigned int compIdx) const
151 { return switchFluxSign_[phaseIdx][compIdx]; }
152
153
154 private:
156 std::array< std::array<bool, numComponents>, numPhases > switchFluxSign_;
157 std::array< std::array<const Stencil*, numComponents>, numPhases > stencil_;
158
160 const PrimaryDataHandle* primaryHandlePtr_;
161 const SecondaryDataHandle* secondaryHandlePtr_;
162 };
163
164public:
165 // state the discretization method this implementation belongs to
167
168 //return the reference system
170 { return referenceSystem; }
171
172 // state the type for the corresponding cache and its filler
173 using Cache = MpfaFicksLawCache;
174
175 // export the diffusion container
177
179 static ComponentFluxVector flux(const Problem& problem,
180 const Element& element,
181 const FVElementGeometry& fvGeometry,
182 const ElementVolumeVariables& elemVolVars,
183 const SubControlVolumeFace& scvf,
184 const int phaseIdx,
185 const ElementFluxVariablesCache& elemFluxVarsCache)
186 {
187 // obtain this scvf's cache
188 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
189
190 ComponentFluxVector componentFlux(0.0);
191 for (int compIdx = 0; compIdx < numComponents; compIdx++)
192 {
193 if constexpr (!FluidSystem::isTracerFluidSystem())
194 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
195 continue;
196
197 // calculate the density at the interface
198 const auto rho = interpolateDensity(elemVolVars, scvf, phaseIdx);
199
200 // compute the flux
201 if (fluxVarsCache.usesSecondaryIv())
202 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
203 fluxVarsCache,
204 fluxVarsCache.diffusionSecondaryDataHandle(),
205 phaseIdx, compIdx);
206 else
207 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
208 fluxVarsCache,
209 fluxVarsCache.diffusionPrimaryDataHandle(),
210 phaseIdx, compIdx);
211 }
212
213 // accumulate the phase component flux
214 for (int compIdx = 0; compIdx < numComponents; compIdx++)
215 if constexpr (!FluidSystem::isTracerFluidSystem())
216 if (compIdx != FluidSystem::getMainComponent(phaseIdx) && BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
217 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
218
219 return componentFlux;
220 }
221
222private:
223 template< class Problem, class FluxVarsCache, class DataHandle >
224 static Scalar computeVolumeFlux(const Problem& problem,
225 const FluxVarsCache& cache,
226 const DataHandle& dataHandle,
227 int phaseIdx, int compIdx)
228 {
229 dataHandle.setPhaseIndex(phaseIdx);
230 dataHandle.setComponentIndex(compIdx);
231
232 const bool switchSign = cache.diffusionSwitchFluxSign(phaseIdx, compIdx);
233
234 const auto localFaceIdx = cache.ivLocalFaceIndex();
235 const auto idxInOutside = cache.indexInOutsideFaces();
236 const auto& xj = dataHandle.uj();
237 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
238 : (!switchSign ? dataHandle.T()[localFaceIdx]
239 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
240 Scalar scvfFlux = tij*xj;
241
242 // switch the sign if necessary
243 if (switchSign)
244 scvfFlux *= -1.0;
245
246 return scvfFlux;
247 }
248
250 static Scalar interpolateDensity(const ElementVolumeVariables& elemVolVars,
251 const SubControlVolumeFace& scvf,
252 const unsigned int phaseIdx)
253 {
254 // use arithmetic mean of the densities around the scvf
255 if (!scvf.boundary())
256 {
257 const Scalar rhoInside = massOrMolarDensity(elemVolVars[scvf.insideScvIdx()], referenceSystem, phaseIdx);
258
259 Scalar rho = rhoInside;
260 for (const auto outsideIdx : scvf.outsideScvIndices())
261 {
262 const Scalar rhoOutside = massOrMolarDensity(elemVolVars[outsideIdx], referenceSystem, phaseIdx);
263 rho += rhoOutside;
264 }
265 return rho/(scvf.outsideScvIndices().size()+1);
266 }
267 else
268 return massOrMolarDensity(elemVolVars[scvf.outsideScvIdx()], referenceSystem, phaseIdx);
269 }
270};
271
272} // end namespace
273
274#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...
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
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: flux/box/fickslaw.hh:43
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/ccmpfa/fickslaw.hh:169
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: flux/ccmpfa/fickslaw.hh:179
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.