24#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_FICKS_LAW_HH
36template<
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
37class FicksLawImplementation;
43template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
49 using Element =
typename GridView::template Codim<0>::Entity;
51 static constexpr int dim = GridView::dimension;
52 static constexpr int dimWorld = GridView::dimensionworld;
55 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
64 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
67 class MpfaFicksLawCacheFiller
72 template<
class FluxVariablesCacheFiller>
73 static void fill(FluxVariablesCache& scvfFluxVarsCache,
74 unsigned int phaseIdx,
unsigned int compIdx,
75 const Problem& problem,
76 const Element& element,
77 const FVElementGeometry& fvGeometry,
78 const ElementVolumeVariables& elemVolVars,
79 const SubControlVolumeFace& scvf,
80 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
83 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
84 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.secondaryInteractionVolume(),
85 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
86 fluxVarsCacheFiller.secondaryIvDataHandle(),
89 scvfFluxVarsCache.updateDiffusion(fluxVarsCacheFiller.primaryInteractionVolume(),
90 fluxVarsCacheFiller.primaryIvLocalFaceData(),
91 fluxVarsCacheFiller.primaryIvDataHandle(),
97 class MpfaFicksLawCache
100 using Stencil =
typename DualGridNodalIndexSet::NodalGridStencilType;
103 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
104 using PrimaryDataHandle =
typename ElementFluxVariablesCache::PrimaryIvDataHandle::DiffusionHandle;
105 using SecondaryDataHandle =
typename ElementFluxVariablesCache::SecondaryIvDataHandle::DiffusionHandle;
108 template<
bool doSecondary = cons
iderSecondaryIVs, std::enable_if_t<doSecondary,
int> = 0 >
109 void setHandlePointer_(
const SecondaryDataHandle& dataHandle)
110 { secondaryHandlePtr_ = &dataHandle; }
113 void setHandlePointer_(
const PrimaryDataHandle& dataHandle)
114 { primaryHandlePtr_ = &dataHandle; }
118 using Filler = MpfaFicksLawCacheFiller;
128 template<
class IV,
class LocalFaceData,
class DataHandle>
129 void updateDiffusion(
const IV& iv,
130 const LocalFaceData& localFaceData,
131 const DataHandle& dataHandle,
132 unsigned int phaseIdx,
unsigned int compIdx)
134 stencil_[phaseIdx][compIdx] = &iv.stencil();
135 switchFluxSign_[phaseIdx][compIdx] = localFaceData.isOutsideFace();
136 setHandlePointer_(dataHandle.diffusionHandle());
140 const Stencil& diffusionStencil(
unsigned int phaseIdx,
unsigned int compIdx)
const
141 {
return *stencil_[phaseIdx][compIdx]; }
144 const PrimaryDataHandle& diffusionPrimaryDataHandle()
const {
return *primaryHandlePtr_; }
145 const SecondaryDataHandle& diffusionSecondaryDataHandle()
const {
return *secondaryHandlePtr_; }
148 bool diffusionSwitchFluxSign(
unsigned int phaseIdx,
unsigned int compIdx)
const
149 {
return switchFluxSign_[phaseIdx][compIdx]; }
154 std::array< std::array<bool, numComponents>, numPhases > switchFluxSign_;
155 std::array< std::array<const Stencil*, numComponents>, numPhases > stencil_;
158 const PrimaryDataHandle* primaryHandlePtr_;
159 const SecondaryDataHandle* secondaryHandlePtr_;
167 {
return referenceSystem; }
172 static ComponentFluxVector
flux(
const Problem& problem,
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const SubControlVolumeFace& scvf,
178 const ElementFluxVariablesCache& elemFluxVarsCache)
181 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
183 ComponentFluxVector componentFlux(0.0);
184 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
186 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
190 const auto effFactor = computeEffectivityFactor(elemVolVars, scvf, phaseIdx);
193 if (effFactor == 0.0)
197 const auto rho = interpolateDensity(elemVolVars, scvf, phaseIdx);
200 if (fluxVarsCache.usesSecondaryIv())
201 componentFlux[compIdx] = rho*effFactor*computeVolumeFlux(problem,
203 fluxVarsCache.diffusionSecondaryDataHandle(),
206 componentFlux[compIdx] = rho*effFactor*computeVolumeFlux(problem,
208 fluxVarsCache.diffusionPrimaryDataHandle(),
213 for(
int compIdx = 0; compIdx < numComponents; compIdx++)
214 if(compIdx != FluidSystem::getMainComponent(phaseIdx) && BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
215 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
217 return componentFlux;
221 template<
class Problem,
class FluxVarsCache,
class DataHandle >
222 static Scalar computeVolumeFlux(
const Problem& problem,
223 const FluxVarsCache& cache,
224 const DataHandle& dataHandle,
225 int phaseIdx,
int compIdx)
227 dataHandle.setPhaseIndex(phaseIdx);
228 dataHandle.setComponentIndex(compIdx);
230 const bool switchSign = cache.diffusionSwitchFluxSign(phaseIdx, compIdx);
232 const auto localFaceIdx = cache.ivLocalFaceIndex();
233 const auto idxInOutside = cache.indexInOutsideFaces();
234 const auto& xj = dataHandle.uj();
235 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
236 : (!switchSign ? dataHandle.T()[localFaceIdx]
237 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
238 Scalar scvfFlux = tij*xj;
248 static Scalar interpolateDensity(
const ElementVolumeVariables& elemVolVars,
249 const SubControlVolumeFace& scvf,
250 const unsigned int phaseIdx)
253 if (!scvf.boundary())
255 const Scalar rhoInside =
massOrMolarDensity(elemVolVars[scvf.insideScvIdx()], referenceSystem, phaseIdx);
257 Scalar rho = rhoInside;
258 for (
const auto outsideIdx : scvf.outsideScvIndices())
260 const Scalar rhoOutside =
massOrMolarDensity(elemVolVars[outsideIdx], referenceSystem, phaseIdx);
263 return rho/(scvf.outsideScvIndices().size()+1);
266 return massOrMolarDensity(elemVolVars[scvf.outsideScvIdx()], referenceSystem, phaseIdx);
273 static Scalar computeEffectivityFactor(
const ElementVolumeVariables& elemVolVars,
274 const SubControlVolumeFace& scvf,
275 const unsigned int phaseIdx)
277 using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
280 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
281 const auto factor = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
282 insideVolVars.saturation(phaseIdx),
285 if (!scvf.boundary())
288 Scalar outsideFactor = 0.0;
289 for (
const auto outsideIdx : scvf.outsideScvIndices())
291 const auto& outsideVolVars = elemVolVars[outsideIdx];
292 outsideFactor += EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
293 outsideVolVars.saturation(phaseIdx),
296 outsideFactor /= scvf.outsideScvIndices().size();
Define some often used mathematical functions.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:68
make the local view function available whenever we use the grid geometry
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: box/fickslaw.hh:38
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: ccmpfa/fickslaw.hh:166
MpfaFicksLawCache Cache
Definition: 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: ccmpfa/fickslaw.hh:172
Declares all properties used in Dumux.