24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
27#include <dune/common/fvector.hh>
42template<
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
43class FicksLawImplementation;
49template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
56 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
59 using GridView =
typename GridGeometry::GridView;
61 using VolumeVariables =
typename ElementVolumeVariables::VolumeVariables;
62 using Element =
typename GridView::template Codim<0>::Entity;
68 static const int dim = GridView::dimension;
69 static const int dimWorld = GridView::dimensionworld;
70 static const int numPhases = ModelTraits::numFluidPhases();
71 static const int numComponents = ModelTraits::numFluidComponents();
73 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
76 class TpfaFicksLawCacheFiller
81 template<
class FluxVariablesCacheFiller>
82 static void fill(FluxVariablesCache& scvfFluxVarsCache,
83 unsigned int phaseIdx,
unsigned int compIdx,
84 const Problem& problem,
85 const Element& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const SubControlVolumeFace& scvf,
89 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
91 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
96 class TpfaFicksLawCache
99 using Filler = TpfaFicksLawCacheFiller;
101 void updateDiffusion(
const Problem& problem,
102 const Element& element,
103 const FVElementGeometry& fvGeometry,
104 const ElementVolumeVariables& elemVolVars,
105 const SubControlVolumeFace &scvf,
106 const unsigned int phaseIdx,
107 const unsigned int compIdx)
109 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
112 const Scalar& diffusionTij(
unsigned int phaseIdx,
unsigned int compIdx)
const
113 {
return tij_[phaseIdx][compIdx]; }
116 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
126 {
return referenceSystem; }
143 static ComponentFluxVector
flux(
const Problem& problem,
144 const Element& element,
145 const FVElementGeometry& fvGeometry,
146 const VolumeVariables& insideVolVars,
147 const VolumeVariables& outsideVolVars,
148 const SubControlVolumeFace& scvf,
150 const ElementFluxVariablesCache& elemFluxVarsCache)
152 if constexpr (isMixedDimensional_)
153 if (scvf.numOutsideScv() != 1)
154 DUNE_THROW(Dune::Exception,
"This flux overload requires scvf.numOutsideScv() == 1");
157 const auto getOutsideX = [&](
const Scalar xInside,
const Scalar tij,
const int phaseIdx,
const int compIdx)
163 const auto getRho = [&](
const int phaseIdx,
const Scalar rhoInside,
const Scalar rhoOutside)
165 return 0.5*(rhoInside + rhoOutside);
168 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
178 static ComponentFluxVector
flux(
const Problem& problem,
179 const Element& element,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars,
182 const SubControlVolumeFace& scvf,
184 const ElementFluxVariablesCache& elemFluxVarsCache)
187 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
188 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
191 const auto getOutsideX = [&](
const Scalar xInside,
const Scalar tij,
const int phaseIdx,
const int compIdx)
193 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
194 if constexpr (isMixedDimensional_)
196 return scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
197 : branchingFacetX_(problem, element, fvGeometry, elemVolVars,
198 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
201 return massOrMoleFractionOutside;
205 const auto getRho = [&](
const int phaseIdx,
const Scalar rhoInside,
const Scalar rhoOutside)
207 if constexpr (isMixedDimensional_)
209 return scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
210 : branchingFacetDensity_(elemVolVars, scvf, phaseIdx, rhoInside);
213 return 0.5*(rhoInside + rhoOutside);
216 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
221 const Element& element,
222 const FVElementGeometry& fvGeometry,
223 const ElementVolumeVariables& elemVolVars,
224 const SubControlVolumeFace& scvf,
225 const int phaseIdx,
const int compIdx)
229 const auto insideScvIdx = scvf.insideScvIdx();
230 const auto& insideScv = fvGeometry.scv(insideScvIdx);
231 const auto& insideVolVars = elemVolVars[insideScvIdx];
232 const auto getDiffCoeff = [&](
const auto& vv)
234 using FluidSystem =
typename ElementVolumeVariables::VolumeVariables::FluidSystem;
235 if constexpr (FluidSystem::isTracerFluidSystem())
236 return vv.effectiveDiffusionCoefficient(0, 0, compIdx);
238 return vv.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
241 const auto insideD = getDiffCoeff(insideVolVars);
247 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
248 tij = Extrusion::area(scvf)*ti;
253 const auto outsideScvIdx = scvf.outsideScvIdx();
254 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
255 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
256 const auto outsideD = getDiffCoeff(outsideVolVars);
259 if constexpr (dim == dimWorld)
266 outsideVolVars.extrusionFactor());
272 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
279 template<
class Outs
ideFractionHelper,
class DensityHelper>
280 static ComponentFluxVector flux_(
const Element& element,
281 const FVElementGeometry& fvGeometry,
282 const VolumeVariables& insideVolVars,
283 const VolumeVariables& outsideVolVars,
284 const ElementFluxVariablesCache& elemFluxVarsCache,
285 const SubControlVolumeFace& scvf,
287 const OutsideFractionHelper& getOutsideX,
288 const DensityHelper& getRho)
290 ComponentFluxVector componentFlux(0.0);
291 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
293 using FluidSystem =
typename VolumeVariables::FluidSystem;
294 if constexpr (!FluidSystem::isTracerFluidSystem())
295 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
299 const Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
302 const Scalar xInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
303 const Scalar xOutside = getOutsideX(xInside, tij, phaseIdx, compIdx);
305 const Scalar rhoInside =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
306 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
308 const Scalar rho = getRho(phaseIdx, rhoInside, rhoOutside);
310 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
311 if constexpr (!FluidSystem::isTracerFluidSystem())
312 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
313 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
316 return componentFlux;
320 static Scalar branchingFacetX_(
const Problem& problem,
321 const Element& element,
322 const FVElementGeometry& fvGeometry,
323 const ElementVolumeVariables& elemVolVars,
324 const ElementFluxVariablesCache& elemFluxVarsCache,
325 const SubControlVolumeFace& scvf,
326 const Scalar insideX,
const Scalar insideTi,
327 const int phaseIdx,
const int compIdx)
329 Scalar sumTi(insideTi);
330 Scalar sumXTi(insideTi*insideX);
332 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
334 const auto outsideScvIdx = scvf.outsideScvIdx(i);
335 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
336 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
337 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
339 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
341 sumXTi += outsideTi*massOrMoleFractionOutside;
344 return sumTi > 0 ? sumXTi/sumTi : 0;
348 static Scalar branchingFacetDensity_(
const ElementVolumeVariables& elemVolVars,
349 const SubControlVolumeFace& scvf,
351 const Scalar insideRho)
353 Scalar rho(insideRho);
354 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
356 const auto outsideScvIdx = scvf.outsideScvIdx(i);
357 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
358 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
361 return rho/(scvf.numOutsideScvs()+1);
364 static constexpr bool isMixedDimensional_ =
static_cast<int>(GridView::dimension) <
static_cast<int>(GridView::dimensionworld);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
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...
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:47
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:66
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 Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fickslaw.hh:51
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/cctpfa/fickslaw.hh:178
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:125
static Scalar calculateTransmissibility(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const int compIdx)
compute diffusive transmissibilities
Definition: flux/cctpfa/fickslaw.hh:220
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:129
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &insideVolVars, const VolumeVariables &outsideVolVars, 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/cctpfa/fickslaw.hh:143
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.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...