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;
64 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
65 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
69 static const int dim = GridView::dimension;
70 static const int dimWorld = GridView::dimensionworld;
71 static const int numPhases = ModelTraits::numFluidPhases();
72 static const int numComponents = ModelTraits::numFluidComponents();
74 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
77 class TpfaFicksLawCacheFiller
82 template<
class FluxVariablesCacheFiller>
83 static void fill(FluxVariablesCache& scvfFluxVarsCache,
84 unsigned int phaseIdx,
unsigned int compIdx,
85 const Problem& problem,
86 const Element& element,
87 const FVElementGeometry& fvGeometry,
88 const ElementVolumeVariables& elemVolVars,
89 const SubControlVolumeFace& scvf,
90 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
92 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
97 class TpfaFicksLawCache
100 using Filler = TpfaFicksLawCacheFiller;
102 void updateDiffusion(
const Problem& problem,
103 const Element& element,
104 const FVElementGeometry& fvGeometry,
105 const ElementVolumeVariables& elemVolVars,
106 const SubControlVolumeFace &scvf,
107 const unsigned int phaseIdx,
108 const unsigned int compIdx)
110 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
113 const Scalar& diffusionTij(
unsigned int phaseIdx,
unsigned int compIdx)
const
114 {
return tij_[phaseIdx][compIdx]; }
117 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
127 {
return referenceSystem; }
144 static ComponentFluxVector
flux(
const Problem& problem,
145 const Element& element,
146 const FVElementGeometry& fvGeometry,
147 const VolumeVariables& insideVolVars,
148 const VolumeVariables& outsideVolVars,
149 const SubControlVolumeFace& scvf,
151 const ElementFluxVariablesCache& elemFluxVarsCache)
153 if constexpr (isMixedDimensional_)
154 if (scvf.numOutsideScv() != 1)
155 DUNE_THROW(Dune::Exception,
"This flux overload requires scvf.numOutsideScv() == 1");
158 const auto getOutsideX = [&](
const Scalar xInside,
const Scalar tij,
const int phaseIdx,
const int compIdx)
164 const auto getRho = [&](
const int phaseIdx,
const Scalar rhoInside,
const Scalar rhoOutside)
166 return 0.5*(rhoInside + rhoOutside);
169 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
179 static ComponentFluxVector
flux(
const Problem& problem,
180 const Element& element,
181 const FVElementGeometry& fvGeometry,
182 const ElementVolumeVariables& elemVolVars,
183 const SubControlVolumeFace& scvf,
185 const ElementFluxVariablesCache& elemFluxVarsCache)
188 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
189 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
192 const auto getOutsideX = [&](
const Scalar xInside,
const Scalar tij,
const int phaseIdx,
const int compIdx)
194 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
195 if constexpr (isMixedDimensional_)
197 return scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
198 : branchingFacetX_(problem, element, fvGeometry, elemVolVars,
199 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
202 return massOrMoleFractionOutside;
206 const auto getRho = [&](
const int phaseIdx,
const Scalar rhoInside,
const Scalar rhoOutside)
208 if constexpr (isMixedDimensional_)
210 return scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
211 : branchingFacetDensity_(elemVolVars, scvf, phaseIdx, rhoInside);
214 return 0.5*(rhoInside + rhoOutside);
217 return flux_(element, fvGeometry, insideVolVars, outsideVolVars, elemFluxVarsCache, scvf, phaseIdx, getOutsideX, getRho);
222 const Element& element,
223 const FVElementGeometry& fvGeometry,
224 const ElementVolumeVariables& elemVolVars,
225 const SubControlVolumeFace& scvf,
226 const int phaseIdx,
const int compIdx)
230 const auto insideScvIdx = scvf.insideScvIdx();
231 const auto& insideScv = fvGeometry.scv(insideScvIdx);
232 const auto& insideVolVars = elemVolVars[insideScvIdx];
233 const auto getDiffCoeff = [&](
const auto& vv)
235 using FluidSystem =
typename ElementVolumeVariables::VolumeVariables::FluidSystem;
236 if constexpr (FluidSystem::isTracerFluidSystem())
237 return vv.effectiveDiffusionCoefficient(0, 0, compIdx);
239 return vv.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
242 const auto insideDiffCoeff = getDiffCoeff(insideVolVars);
248 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
249 tij = Extrusion::area(fvGeometry, scvf)*ti;
254 const auto outsideScvIdx = scvf.outsideScvIdx();
255 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
256 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
257 const auto outsideDiffCoeff = getDiffCoeff(outsideVolVars);
260 if constexpr (dim == dimWorld)
267 outsideVolVars.extrusionFactor());
273 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
280 template<
class Outs
ideFractionHelper,
class DensityHelper>
281 static ComponentFluxVector flux_(
const Element& element,
282 const FVElementGeometry& fvGeometry,
283 const VolumeVariables& insideVolVars,
284 const VolumeVariables& outsideVolVars,
285 const ElementFluxVariablesCache& elemFluxVarsCache,
286 const SubControlVolumeFace& scvf,
288 const OutsideFractionHelper& getOutsideX,
289 const DensityHelper& getRho)
291 ComponentFluxVector componentFlux(0.0);
292 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
294 using FluidSystem =
typename VolumeVariables::FluidSystem;
295 if constexpr (!FluidSystem::isTracerFluidSystem())
296 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
300 const Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
303 const Scalar xInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
304 const Scalar xOutside = getOutsideX(xInside, tij, phaseIdx, compIdx);
306 const Scalar rhoInside =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
307 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
309 const Scalar rho = getRho(phaseIdx, rhoInside, rhoOutside);
311 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
312 if constexpr (!FluidSystem::isTracerFluidSystem())
313 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
314 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
317 return componentFlux;
321 static Scalar branchingFacetX_(
const Problem& problem,
322 const Element& element,
323 const FVElementGeometry& fvGeometry,
324 const ElementVolumeVariables& elemVolVars,
325 const ElementFluxVariablesCache& elemFluxVarsCache,
326 const SubControlVolumeFace& scvf,
327 const Scalar insideX,
const Scalar insideTi,
328 const int phaseIdx,
const int compIdx)
330 Scalar sumTi(insideTi);
331 Scalar sumXTi(insideTi*insideX);
333 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
335 const auto outsideScvIdx = scvf.outsideScvIdx(i);
336 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
337 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
338 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
340 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
342 sumXTi += outsideTi*massOrMoleFractionOutside;
345 return sumTi > 0 ? sumXTi/sumTi : 0;
349 static Scalar branchingFacetDensity_(
const ElementVolumeVariables& elemVolVars,
350 const SubControlVolumeFace& scvf,
352 const Scalar insideRho)
354 Scalar rho(insideRho);
355 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
357 const auto outsideScvIdx = scvf.outsideScvIdx(i);
358 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
359 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
362 return rho/(scvf.numOutsideScvs()+1);
365 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:48
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
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
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:179
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:126
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:221
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:130
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:144
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...