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, DiscretizationMethod discMethod, 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_;
170 {
return referenceSystem; }
184 static ComponentFluxVector
flux(
const Problem& problem,
185 const Element& element,
186 const FVElementGeometry& fvGeometry,
187 const ElementVolumeVariables& elemVolVars,
188 const SubControlVolumeFace& scvf,
190 const ElementFluxVariablesCache& elemFluxVarsCache)
193 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
195 ComponentFluxVector componentFlux(0.0);
196 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
198 if constexpr (!FluidSystem::isTracerFluidSystem())
199 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
203 const auto rho = interpolateDensity(elemVolVars, scvf, phaseIdx);
206 if (fluxVarsCache.usesSecondaryIv())
207 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
209 fluxVarsCache.diffusionSecondaryDataHandle(),
212 componentFlux[compIdx] = rho*computeVolumeFlux(problem,
214 fluxVarsCache.diffusionPrimaryDataHandle(),
219 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
220 if constexpr (!FluidSystem::isTracerFluidSystem())
221 if (compIdx != FluidSystem::getMainComponent(phaseIdx) && BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
222 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
224 return componentFlux;
228 template<
class Problem,
class FluxVarsCache,
class DataHandle >
229 static Scalar computeVolumeFlux(
const Problem& problem,
230 const FluxVarsCache& cache,
231 const DataHandle& dataHandle,
232 int phaseIdx,
int compIdx)
234 dataHandle.setPhaseIndex(phaseIdx);
235 dataHandle.setComponentIndex(compIdx);
237 const bool switchSign = cache.diffusionSwitchFluxSign(phaseIdx, compIdx);
239 const auto localFaceIdx = cache.ivLocalFaceIndex();
240 const auto idxInOutside = cache.indexInOutsideFaces();
241 const auto& xj = dataHandle.uj();
242 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
243 : (!switchSign ? dataHandle.T()[localFaceIdx]
244 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
245 Scalar scvfFlux = tij*xj;
255 static Scalar interpolateDensity(
const ElementVolumeVariables& elemVolVars,
256 const SubControlVolumeFace& scvf,
257 const unsigned int phaseIdx)
260 if (!scvf.boundary())
262 const Scalar rhoInside =
massOrMolarDensity(elemVolVars[scvf.insideScvIdx()], referenceSystem, phaseIdx);
264 Scalar rho = rhoInside;
265 for (
const auto outsideIdx : scvf.outsideScvIndices())
267 const Scalar rhoOutside =
massOrMolarDensity(elemVolVars[outsideIdx], referenceSystem, phaseIdx);
270 return rho/(scvf.outsideScvIndices().size()+1);
273 return massOrMolarDensity(elemVolVars[scvf.outsideScvIdx()], referenceSystem, phaseIdx);
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...
The available discretization methods in Dumux.
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
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:43
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: flux/ccmpfa/fickslaw.hh:169
MpfaFicksLawCache Cache
Definition: flux/ccmpfa/fickslaw.hh:173
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:184
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.