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;
66 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
69 class MpfaFicksLawCacheFiller
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)
85 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
86 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.secondaryInteractionVolume(),
87 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
88 fluxVarsCacheFiller.secondaryIvDataHandle(),
91 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.primaryInteractionVolume(),
92 fluxVarsCacheFiller.primaryIvLocalFaceData(),
93 fluxVarsCacheFiller.primaryIvDataHandle(),
99 class MpfaFicksLawCache
102 using Stencil =
typename DualGridNodalIndexSet::NodalGridStencilType;
105 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
106 using PrimaryDataHandle =
typename ElementFluxVariablesCache::PrimaryIvDataHandle::DiffusionHandle;
107 using SecondaryDataHandle =
typename ElementFluxVariablesCache::SecondaryIvDataHandle::DiffusionHandle;
110 template<
bool doSecondary = cons
iderSecondaryIVs, std::enable_if_t<doSecondary,
int> = 0 >
111 void setHandlePointer_(
const SecondaryDataHandle& dataHandle)
112 { secondaryHandlePtr_ = &dataHandle; }
115 void setHandlePointer_(
const PrimaryDataHandle& dataHandle)
116 { primaryHandlePtr_ = &dataHandle; }
120 using Filler = MpfaFicksLawCacheFiller;
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)
136 stencil_[phaseIdx][compIdx] = &iv.stencil();
137 switchFluxSign_[phaseIdx][compIdx] = localFaceData.isOutsideFace();
138 setHandlePointer_(dataHandle.diffusionHandle());
142 const Stencil& diffusionStencil(
unsigned int phaseIdx,
unsigned int compIdx)
const
143 {
return *stencil_[phaseIdx][compIdx]; }
146 const PrimaryDataHandle& diffusionPrimaryDataHandle()
const {
return *primaryHandlePtr_; }
147 const SecondaryDataHandle& diffusionSecondaryDataHandle()
const {
return *secondaryHandlePtr_; }
150 bool diffusionSwitchFluxSign(
unsigned int phaseIdx,
unsigned int compIdx)
const
151 {
return switchFluxSign_[phaseIdx][compIdx]; }
156 std::array< std::array<bool, numComponents>, numPhases > switchFluxSign_;
157 std::array< std::array<const Stencil*, numComponents>, numPhases > stencil_;
160 const PrimaryDataHandle* primaryHandlePtr_;
161 const SecondaryDataHandle* secondaryHandlePtr_;
171 {
return referenceSystem; }
185 static ComponentFluxVector
flux(
const Problem& problem,
186 const Element& element,
187 const FVElementGeometry& fvGeometry,
188 const ElementVolumeVariables& elemVolVars,
189 const SubControlVolumeFace& scvf,
191 const ElementFluxVariablesCache& elemFluxVarsCache)
194 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
196 ComponentFluxVector componentFlux(0.0);
197 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
199 if constexpr (!FluidSystem::isTracerFluidSystem())
200 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
204 const auto rho = interpolateDensity(elemVolVars, scvf, phaseIdx);
207 if (fluxVarsCache.usesSecondaryIv())
208 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
210 fluxVarsCache.diffusionSecondaryDataHandle(),
213 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
215 fluxVarsCache.diffusionPrimaryDataHandle(),
220 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
221 if constexpr (!FluidSystem::isTracerFluidSystem())
222 if (compIdx != FluidSystem::getMainComponent(phaseIdx) && BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
223 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
225 return componentFlux;
229 template<
class Problem,
class FluxVarsCache,
class DataHandle >
230 static Scalar computeVolumeFlux(
const Problem& problem,
231 const FluxVarsCache& cache,
232 const DataHandle& dataHandle,
233 int phaseIdx,
int compIdx)
235 dataHandle.setPhaseIndex(phaseIdx);
236 dataHandle.setComponentIndex(compIdx);
238 const bool switchSign = cache.diffusionSwitchFluxSign(phaseIdx, compIdx);
240 const auto localFaceIdx = cache.ivLocalFaceIndex();
241 const auto idxInOutside = cache.indexInOutsideFaces();
242 const auto& xj = dataHandle.uj();
243 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
244 : (!switchSign ? dataHandle.T()[localFaceIdx]
245 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
246 Scalar scvfFlux = tij*xj;
256 static Scalar interpolateDensity(
const ElementVolumeVariables& elemVolVars,
257 const SubControlVolumeFace& scvf,
258 const unsigned int phaseIdx)
261 if (!scvf.boundary())
263 const Scalar rhoInside =
massOrMolarDensity(elemVolVars[scvf.insideScvIdx()], referenceSystem, phaseIdx);
265 Scalar rho = rhoInside;
266 for (
const auto outsideIdx : scvf.outsideScvIndices())
268 const Scalar rhoOutside =
massOrMolarDensity(elemVolVars[outsideIdx], referenceSystem, phaseIdx);
271 return rho/(scvf.outsideScvIndices().size()+1);
274 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
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
forward declaration of the method-specific implemetation
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:185
MpfaFicksLawCache Cache
Definition: flux/ccmpfa/fickslaw.hh:174
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/ccmpfa/fickslaw.hh:170
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.