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 SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
59 using Element =
typename GridView::template Codim<0>::Entity;
66 static const int dim = GridView::dimension;
67 static const int dimWorld = GridView::dimensionworld;
68 static const int numPhases = ModelTraits::numFluidPhases();
69 static const int numComponents = ModelTraits::numFluidComponents();
71 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
74 class TpfaFicksLawCacheFiller
79 template<
class FluxVariablesCacheFiller>
80 static void fill(FluxVariablesCache& scvfFluxVarsCache,
81 unsigned int phaseIdx,
unsigned int compIdx,
82 const Problem& problem,
83 const Element& element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const SubControlVolumeFace& scvf,
87 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
89 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
94 class TpfaFicksLawCache
97 using Filler = TpfaFicksLawCacheFiller;
99 void updateDiffusion(
const Problem& problem,
100 const Element& element,
101 const FVElementGeometry& fvGeometry,
102 const ElementVolumeVariables& elemVolVars,
103 const SubControlVolumeFace &scvf,
104 const unsigned int phaseIdx,
105 const unsigned int compIdx)
107 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
110 const Scalar& diffusionTij(
unsigned int phaseIdx,
unsigned int compIdx)
const
111 {
return tij_[phaseIdx][compIdx]; }
114 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
122 {
return referenceSystem; }
130 static ComponentFluxVector
flux(
const Problem& problem,
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& elemVolVars,
134 const SubControlVolumeFace& scvf,
136 const ElementFluxVariablesCache& elemFluxVarsCache)
138 ComponentFluxVector componentFlux(0.0);
139 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
141 if constexpr (!FluidSystem::isTracerFluidSystem())
142 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
146 Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
149 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
150 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
153 const Scalar xInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
154 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
155 const Scalar xOutside = scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
156 : branchingFacetX(problem, element, fvGeometry, elemVolVars,
157 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
159 const Scalar rhoInside =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
160 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
162 const Scalar rho = scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
163 : branchingFacetDensity(elemVolVars, scvf, phaseIdx, rhoInside);
165 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
166 if constexpr (!FluidSystem::isTracerFluidSystem())
167 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
168 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
171 return componentFlux;
176 const Element& element,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const SubControlVolumeFace& scvf,
180 const int phaseIdx,
const int compIdx)
184 const auto insideScvIdx = scvf.insideScvIdx();
185 const auto& insideScv = fvGeometry.scv(insideScvIdx);
186 const auto& insideVolVars = elemVolVars[insideScvIdx];
187 const auto getDiffCoeff = [&](
const auto& vv)
190 if constexpr (FluidSystem::isTracerFluidSystem())
191 return Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(vv, 0, 0, compIdx);
193 return Deprecated::template effectiveDiffusionCoefficient<EffDiffModel>(vv, phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
196 const auto insideD = getDiffCoeff(insideVolVars);
202 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
203 tij = scvf.area()*ti;
208 const auto outsideScvIdx = scvf.outsideScvIdx();
209 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
210 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
211 const auto outsideD = getDiffCoeff(outsideVolVars);
221 outsideVolVars.extrusionFactor());
227 tij = scvf.area()*(ti * tj)/(ti + tj);
235 static Scalar branchingFacetX(
const Problem& problem,
236 const Element& element,
237 const FVElementGeometry& fvGeometry,
238 const ElementVolumeVariables& elemVolVars,
239 const ElementFluxVariablesCache& elemFluxVarsCache,
240 const SubControlVolumeFace& scvf,
241 const Scalar insideX,
const Scalar insideTi,
242 const int phaseIdx,
const int compIdx)
244 Scalar sumTi(insideTi);
245 Scalar sumXTi(insideTi*insideX);
247 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
249 const auto outsideScvIdx = scvf.outsideScvIdx(i);
250 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
251 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
252 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
254 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
256 sumXTi += outsideTi*massOrMoleFractionOutside;
259 return sumTi > 0 ? sumXTi/sumTi : 0;
263 static Scalar branchingFacetDensity(
const ElementVolumeVariables& elemVolVars,
264 const SubControlVolumeFace& scvf,
266 const Scalar insideRho)
268 Scalar rho(insideRho);
269 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
271 const auto outsideScvIdx = scvf.outsideScvIdx(i);
272 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
273 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
276 return rho/(scvf.numOutsideScvs()+1);
Container storing the diffusion coefficients required by Fick's law. Uses the minimal possible contai...
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
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 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:175
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: flux/cctpfa/fickslaw.hh:121
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: flux/cctpfa/fickslaw.hh:125
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
return diffusive fluxes for all components in a phase
Definition: flux/cctpfa/fickslaw.hh:130
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...