version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FICKS_LAW_HH
13#define DUMUX_DISCRETIZATION_CC_MPFA_FICKS_LAW_HH
14
15#include <dune/common/fvector.hh>
18
21
22namespace Dumux {
23
25template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
26class FicksLawImplementation;
27
32template <class TypeTag, ReferenceSystemFormulation referenceSystem>
33class FicksLawImplementation<TypeTag, DiscretizationMethods::CCMpfa, referenceSystem>
34{
38 using Element = typename GridView::template Codim<0>::Entity;
39
40 static constexpr int dim = GridView::dimension;
41 static constexpr int dimWorld = GridView::dimensionworld;
42
44 using FVElementGeometry = typename GridGeometry::LocalView;
46 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
47 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
49 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
50 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
52
53 static constexpr int numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
55 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
56
58 class MpfaFicksLawCacheFiller
59 {
60 public:
63 template<class FluxVariablesCacheFiller>
64 static void fill(FluxVariablesCache& scvfFluxVarsCache,
65 unsigned int phaseIdx, unsigned int compIdx,
66 const Problem& problem,
67 const Element& element,
68 const FVElementGeometry& fvGeometry,
69 const ElementVolumeVariables& elemVolVars,
70 const SubControlVolumeFace& scvf,
71 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
72 {
73 // get interaction volume related data from the filler class & update the cache
74 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
75 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.secondaryInteractionVolume(),
76 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
77 fluxVarsCacheFiller.secondaryIvDataHandle(),
78 phaseIdx, compIdx);
79 else
80 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.primaryInteractionVolume(),
81 fluxVarsCacheFiller.primaryIvLocalFaceData(),
82 fluxVarsCacheFiller.primaryIvDataHandle(),
83 phaseIdx, compIdx);
84 }
85 };
86
88 class MpfaFicksLawCache
89 {
91 using Stencil = typename DualGridNodalIndexSet::NodalGridStencilType;
92
94 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
95 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle::DiffusionHandle;
96 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle::DiffusionHandle;
97
99 template< bool doSecondary = considerSecondaryIVs, std::enable_if_t<doSecondary, int> = 0 >
100 void setHandlePointer_(const SecondaryDataHandle& dataHandle)
101 { secondaryHandlePtr_ = &dataHandle; }
102
104 void setHandlePointer_(const PrimaryDataHandle& dataHandle)
105 { primaryHandlePtr_ = &dataHandle; }
106
107 public:
108 // export filler type
109 using Filler = MpfaFicksLawCacheFiller;
110
119 template<class IV, class LocalFaceData, class DataHandle>
120 void updateDiffusion(const IV& iv,
121 const LocalFaceData& localFaceData,
122 const DataHandle& dataHandle,
123 unsigned int phaseIdx, unsigned int compIdx)
124 {
125 stencil_[phaseIdx][compIdx] = &iv.stencil();
126 switchFluxSign_[phaseIdx][compIdx] = localFaceData.isOutsideFace();
127 setHandlePointer_(dataHandle.diffusionHandle());
128 }
129
131 const Stencil& diffusionStencil(unsigned int phaseIdx, unsigned int compIdx) const
132 { return *stencil_[phaseIdx][compIdx]; }
133
135 const PrimaryDataHandle& diffusionPrimaryDataHandle() const { return *primaryHandlePtr_; }
136 const SecondaryDataHandle& diffusionSecondaryDataHandle() const { return *secondaryHandlePtr_; }
137
139 bool diffusionSwitchFluxSign(unsigned int phaseIdx, unsigned int compIdx) const
140 { return switchFluxSign_[phaseIdx][compIdx]; }
141
142
143 private:
145 std::array< std::array<bool, numComponents>, numPhases > switchFluxSign_;
146 std::array< std::array<const Stencil*, numComponents>, numPhases > stencil_;
147
149 const PrimaryDataHandle* primaryHandlePtr_;
150 const SecondaryDataHandle* secondaryHandlePtr_;
151 };
152
153public:
155 // state the discretization method this implementation belongs to
156 static constexpr DiscretizationMethod discMethod{};
157
158 //return the reference system
160 { return referenceSystem; }
161
162 // state the type for the corresponding cache and its filler
163 using Cache = MpfaFicksLawCache;
164
165 // export the diffusion container
167
174 static ComponentFluxVector flux(const Problem& problem,
175 const Element& element,
176 const FVElementGeometry& fvGeometry,
177 const ElementVolumeVariables& elemVolVars,
178 const SubControlVolumeFace& scvf,
179 const int phaseIdx,
180 const ElementFluxVariablesCache& elemFluxVarsCache)
181 {
182 // obtain this scvf's cache
183 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
184
185 ComponentFluxVector componentFlux(0.0);
186 for (int compIdx = 0; compIdx < numComponents; compIdx++)
187 {
188 if constexpr (!FluidSystem::isTracerFluidSystem())
189 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
190 continue;
191
192 // calculate the density at the interface
193 const auto rho = interpolateDensity(elemVolVars, scvf, phaseIdx);
194
195 // compute the flux
196 if (fluxVarsCache.usesSecondaryIv())
197 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
198 fluxVarsCache,
199 fluxVarsCache.diffusionSecondaryDataHandle(),
200 phaseIdx, compIdx);
201 else
202 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
203 fluxVarsCache,
204 fluxVarsCache.diffusionPrimaryDataHandle(),
205 phaseIdx, compIdx);
206 }
207
208 // accumulate the phase component flux
209 for (int compIdx = 0; compIdx < numComponents; compIdx++)
210 if constexpr (!FluidSystem::isTracerFluidSystem())
211 if (compIdx != FluidSystem::getMainComponent(phaseIdx) && BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
212 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
213
214 return componentFlux;
215 }
216
217private:
218 template< class Problem, class FluxVarsCache, class DataHandle >
219 static Scalar computeVolumeFlux(const Problem& problem,
220 const FluxVarsCache& cache,
221 const DataHandle& dataHandle,
222 int phaseIdx, int compIdx)
223 {
224 dataHandle.setPhaseIndex(phaseIdx);
225 dataHandle.setComponentIndex(compIdx);
226
227 const bool switchSign = cache.diffusionSwitchFluxSign(phaseIdx, compIdx);
228
229 const auto localFaceIdx = cache.ivLocalFaceIndex();
230 const auto idxInOutside = cache.indexInOutsideFaces();
231 const auto& xj = dataHandle.uj();
232 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
233 : (!switchSign ? dataHandle.T()[localFaceIdx]
234 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
235 Scalar scvfFlux = tij*xj;
236
237 // switch the sign if necessary
238 if (switchSign)
239 scvfFlux *= -1.0;
240
241 return scvfFlux;
242 }
243
245 static Scalar interpolateDensity(const ElementVolumeVariables& elemVolVars,
246 const SubControlVolumeFace& scvf,
247 const unsigned int phaseIdx)
248 {
249 // use arithmetic mean of the densities around the scvf
250 if (!scvf.boundary())
251 {
252 const Scalar rhoInside = massOrMolarDensity(elemVolVars[scvf.insideScvIdx()], referenceSystem, phaseIdx);
253
254 Scalar rho = rhoInside;
255 for (const auto outsideIdx : scvf.outsideScvIndices())
256 {
257 const Scalar rhoOutside = massOrMolarDensity(elemVolVars[outsideIdx], referenceSystem, phaseIdx);
258 rho += rhoOutside;
259 }
260 return rho/(scvf.outsideScvIndices().size()+1);
261 }
262 else
263 return massOrMolarDensity(elemVolVars[scvf.outsideScvIdx()], referenceSystem, phaseIdx);
264 }
265};
266
267} // end namespace
268
269#endif
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
Definition: fickiandiffusioncoefficients.hh:32
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:174
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/ccmpfa/fickslaw.hh:159
forward declaration of the method-specific implementation
Definition: flux/box/fickslaw.hh:32
Defines all properties used in Dumux.
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:43
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
Definition: adapt.hh:17
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Definition: method.hh:33