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, DiscretizationMethod discMethod, 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 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_;
124 {
return referenceSystem; }
137 static ComponentFluxVector
flux(
const Problem& problem,
138 const Element& element,
139 const FVElementGeometry& fvGeometry,
140 const ElementVolumeVariables& elemVolVars,
141 const SubControlVolumeFace& scvf,
143 const ElementFluxVariablesCache& elemFluxVarsCache)
145 ComponentFluxVector componentFlux(0.0);
146 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
148 if constexpr (!FluidSystem::isTracerFluidSystem())
149 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
153 Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
156 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
157 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
160 const Scalar xInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
161 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
162 const Scalar xOutside = scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
163 : branchingFacetX(problem, element, fvGeometry, elemVolVars,
164 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
166 const Scalar rhoInside =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
167 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
169 const Scalar rho = scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
170 : branchingFacetDensity(elemVolVars, scvf, phaseIdx, rhoInside);
172 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
173 if constexpr (!FluidSystem::isTracerFluidSystem())
174 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
175 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
178 return componentFlux;
183 const Element& element,
184 const FVElementGeometry& fvGeometry,
185 const ElementVolumeVariables& elemVolVars,
186 const SubControlVolumeFace& scvf,
187 const int phaseIdx,
const int compIdx)
191 const auto insideScvIdx = scvf.insideScvIdx();
192 const auto& insideScv = fvGeometry.scv(insideScvIdx);
193 const auto& insideVolVars = elemVolVars[insideScvIdx];
194 const auto getDiffCoeff = [&](
const auto& vv)
196 if constexpr (FluidSystem::isTracerFluidSystem())
197 return vv.effectiveDiffusionCoefficient(0, 0, compIdx);
199 return vv.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
202 const auto insideD = getDiffCoeff(insideVolVars);
208 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
209 tij = Extrusion::area(scvf)*ti;
214 const auto outsideScvIdx = scvf.outsideScvIdx();
215 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
216 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
217 const auto outsideD = getDiffCoeff(outsideVolVars);
227 outsideVolVars.extrusionFactor());
233 tij = Extrusion::area(scvf)*(ti * tj)/(ti + tj);
241 static Scalar branchingFacetX(
const Problem& problem,
242 const Element& element,
243 const FVElementGeometry& fvGeometry,
244 const ElementVolumeVariables& elemVolVars,
245 const ElementFluxVariablesCache& elemFluxVarsCache,
246 const SubControlVolumeFace& scvf,
247 const Scalar insideX,
const Scalar insideTi,
248 const int phaseIdx,
const int compIdx)
250 Scalar sumTi(insideTi);
251 Scalar sumXTi(insideTi*insideX);
253 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
255 const auto outsideScvIdx = scvf.outsideScvIdx(i);
256 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
257 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
258 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
260 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
262 sumXTi += outsideTi*massOrMoleFractionOutside;
265 return sumTi > 0 ? sumXTi/sumTi : 0;
269 static Scalar branchingFacetDensity(
const ElementVolumeVariables& elemVolVars,
270 const SubControlVolumeFace& scvf,
272 const Scalar insideRho)
274 Scalar rho(insideRho);
275 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
277 const auto outsideScvIdx = scvf.outsideScvIdx(i);
278 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
279 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
282 return rho/(scvf.numOutsideScvs()+1);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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 (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declaration of the method-specific implemetation
Definition: flux/box/fickslaw.hh:43
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: flux/cctpfa/fickslaw.hh:51
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:182
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:123
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:127
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:137
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...