24#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_FICKS_LAW_HH
27#include <dune/common/fvector.hh>
37template<
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
38class FicksLawImplementation;
44template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
50 using Element =
typename GridView::template Codim<0>::Entity;
52 static constexpr int dim = GridView::dimension;
53 static constexpr int dimWorld = GridView::dimensionworld;
56 using FVElementGeometry =
typename GridGeometry::LocalView;
58 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
61 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
62 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
67 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
70 class MpfaFicksLawCacheFiller
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)
86 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
87 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.secondaryInteractionVolume(),
88 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
89 fluxVarsCacheFiller.secondaryIvDataHandle(),
92 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.primaryInteractionVolume(),
93 fluxVarsCacheFiller.primaryIvLocalFaceData(),
94 fluxVarsCacheFiller.primaryIvDataHandle(),
100 class MpfaFicksLawCache
103 using Stencil =
typename DualGridNodalIndexSet::NodalGridStencilType;
106 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
107 using PrimaryDataHandle =
typename ElementFluxVariablesCache::PrimaryIvDataHandle::DiffusionHandle;
108 using SecondaryDataHandle =
typename ElementFluxVariablesCache::SecondaryIvDataHandle::DiffusionHandle;
111 template<
bool doSecondary = cons
iderSecondaryIVs, std::enable_if_t<doSecondary,
int> = 0 >
112 void setHandlePointer_(
const SecondaryDataHandle& dataHandle)
113 { secondaryHandlePtr_ = &dataHandle; }
116 void setHandlePointer_(
const PrimaryDataHandle& dataHandle)
117 { primaryHandlePtr_ = &dataHandle; }
121 using Filler = MpfaFicksLawCacheFiller;
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)
137 stencil_[phaseIdx][compIdx] = &iv.stencil();
138 switchFluxSign_[phaseIdx][compIdx] = localFaceData.isOutsideFace();
139 setHandlePointer_(dataHandle.diffusionHandle());
143 const Stencil& diffusionStencil(
unsigned int phaseIdx,
unsigned int compIdx)
const
144 {
return *stencil_[phaseIdx][compIdx]; }
147 const PrimaryDataHandle& diffusionPrimaryDataHandle()
const {
return *primaryHandlePtr_; }
148 const SecondaryDataHandle& diffusionSecondaryDataHandle()
const {
return *secondaryHandlePtr_; }
151 bool diffusionSwitchFluxSign(
unsigned int phaseIdx,
unsigned int compIdx)
const
152 {
return switchFluxSign_[phaseIdx][compIdx]; }
157 std::array< std::array<bool, numComponents>, numPhases > switchFluxSign_;
158 std::array< std::array<const Stencil*, numComponents>, numPhases > stencil_;
161 const PrimaryDataHandle* primaryHandlePtr_;
162 const SecondaryDataHandle* secondaryHandlePtr_;
172 {
return referenceSystem; }
186 static ComponentFluxVector
flux(
const Problem& problem,
187 const Element& element,
188 const FVElementGeometry& fvGeometry,
189 const ElementVolumeVariables& elemVolVars,
190 const SubControlVolumeFace& scvf,
192 const ElementFluxVariablesCache& elemFluxVarsCache)
195 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
197 ComponentFluxVector componentFlux(0.0);
198 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
200 if constexpr (!FluidSystem::isTracerFluidSystem())
201 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
205 const auto rho = interpolateDensity(elemVolVars, scvf, phaseIdx);
208 if (fluxVarsCache.usesSecondaryIv())
209 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
211 fluxVarsCache.diffusionSecondaryDataHandle(),
214 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
216 fluxVarsCache.diffusionPrimaryDataHandle(),
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];
226 return componentFlux;
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)
236 dataHandle.setPhaseIndex(phaseIdx);
237 dataHandle.setComponentIndex(compIdx);
239 const bool switchSign = cache.diffusionSwitchFluxSign(phaseIdx, compIdx);
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;
257 static Scalar interpolateDensity(
const ElementVolumeVariables& elemVolVars,
258 const SubControlVolumeFace& scvf,
259 const unsigned int phaseIdx)
262 if (!scvf.boundary())
264 const Scalar rhoInside =
massOrMolarDensity(elemVolVars[scvf.insideScvIdx()], referenceSystem, phaseIdx);
266 Scalar rho = rhoInside;
267 for (
const auto outsideIdx : scvf.outsideScvIndices())
269 const Scalar rhoOutside =
massOrMolarDensity(elemVolVars[outsideIdx], referenceSystem, phaseIdx);
272 return rho/(scvf.outsideScvIndices().size()+1);
275 return massOrMolarDensity(elemVolVars[scvf.outsideScvIdx()], referenceSystem, phaseIdx);
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
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
MpfaFicksLawCache Cache
Definition: flux/ccmpfa/fickslaw.hh:175
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.