24#ifndef DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
25#define DUMUX_DISCRETIZATION_CC_TPFA_FICKS_LAW_HH
38template<
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
39class FicksLawImplementation;
45template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
52 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
53 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
57 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 DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
72 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
75 class TpfaFicksLawCacheFiller
80 template<
class FluxVariablesCacheFiller>
81 static void fill(FluxVariablesCache& scvfFluxVarsCache,
82 unsigned int phaseIdx,
unsigned int compIdx,
83 const Problem& problem,
84 const Element& element,
85 const FVElementGeometry& fvGeometry,
86 const ElementVolumeVariables& elemVolVars,
87 const SubControlVolumeFace& scvf,
88 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
90 scvfFluxVarsCache.updateDiffusion(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
95 class TpfaFicksLawCache
98 using Filler = TpfaFicksLawCacheFiller;
100 void updateDiffusion(
const Problem& problem,
101 const Element& element,
102 const FVElementGeometry& fvGeometry,
103 const ElementVolumeVariables& elemVolVars,
104 const SubControlVolumeFace &scvf,
105 const unsigned int phaseIdx,
106 const unsigned int compIdx)
108 tij_[phaseIdx][compIdx] = Implementation::calculateTransmissibility(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, compIdx);
111 const Scalar& diffusionTij(
unsigned int phaseIdx,
unsigned int compIdx)
const
112 {
return tij_[phaseIdx][compIdx]; }
115 std::array< std::array<Scalar, numComponents>, numPhases> tij_;
123 {
return referenceSystem; }
129 static ComponentFluxVector
flux(
const Problem& problem,
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
135 const ElementFluxVariablesCache& elemFluxVarsCache)
137 ComponentFluxVector componentFlux(0.0);
138 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
140 if(compIdx == FluidSystem::getMainComponent(phaseIdx))
144 Scalar tij = elemFluxVarsCache[scvf].diffusionTij(phaseIdx, compIdx);
147 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
148 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
151 const Scalar xInside =
massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
152 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
153 const Scalar xOutside = scvf.numOutsideScvs() == 1 ? massOrMoleFractionOutside
154 : branchingFacetX(problem, element, fvGeometry, elemVolVars,
155 elemFluxVarsCache, scvf, xInside, tij, phaseIdx, compIdx);
157 const Scalar rhoInside =
massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
158 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
160 const Scalar rho = scvf.numOutsideScvs() == 1 ? 0.5*(rhoInside + rhoOutside)
161 : branchingFacetDensity(elemVolVars, scvf, phaseIdx, rhoInside);
163 componentFlux[compIdx] = rho*tij*(xInside - xOutside);
164 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
165 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
168 return componentFlux;
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const SubControlVolumeFace& scvf,
177 const int phaseIdx,
const int compIdx)
181 const auto insideScvIdx = scvf.insideScvIdx();
182 const auto& insideScv = fvGeometry.scv(insideScvIdx);
183 const auto& insideVolVars = elemVolVars[insideScvIdx];
186 const auto insideD = EffDiffModel::effectiveDiffusivity(insideVolVars.porosity(),
187 insideVolVars.saturation(phaseIdx),
188 insideVolVars.diffusionCoefficient(phaseIdx, compIdx));
192 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
194 tij = scvf.area()*ti;
199 const auto outsideScvIdx = scvf.outsideScvIdx();
200 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
201 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
203 const auto outsideD = EffDiffModel::effectiveDiffusivity(outsideVolVars.porosity(),
204 outsideVolVars.saturation(phaseIdx),
205 outsideVolVars.diffusionCoefficient(phaseIdx, compIdx));
215 outsideVolVars.extrusionFactor());
221 tij = scvf.area()*(ti * tj)/(ti + tj);
229 static Scalar branchingFacetX(
const Problem& problem,
230 const Element& element,
231 const FVElementGeometry& fvGeometry,
232 const ElementVolumeVariables& elemVolVars,
233 const ElementFluxVariablesCache& elemFluxVarsCache,
234 const SubControlVolumeFace& scvf,
235 const Scalar insideX,
const Scalar insideTi,
236 const int phaseIdx,
const int compIdx)
238 Scalar sumTi(insideTi);
239 Scalar sumXTi(insideTi*insideX);
241 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
243 const auto outsideScvIdx = scvf.outsideScvIdx(i);
244 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
245 const Scalar massOrMoleFractionOutside =
massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
246 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
248 const Scalar outsideTi = elemFluxVarsCache[flippedScvf].diffusionTij(phaseIdx, compIdx);
250 sumXTi += outsideTi*massOrMoleFractionOutside;
253 return sumTi > 0 ? sumXTi/sumTi : 0;
257 static Scalar branchingFacetDensity(
const ElementVolumeVariables& elemVolVars,
258 const SubControlVolumeFace& scvf,
260 const Scalar insideRho)
262 Scalar rho(insideRho);
263 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
265 const auto outsideScvIdx = scvf.outsideScvIdx(i);
266 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
267 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
270 return rho/(scvf.numOutsideScvs()+1);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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
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
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
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
forward declaration of the method-specific implemetation
Definition: box/fickslaw.hh:38
Fick's law for cell-centered finite volume schemes with two-point flux approximation.
Definition: cctpfa/fickslaw.hh:47
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: cctpfa/fickslaw.hh:172
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Return the reference system.
Definition: cctpfa/fickslaw.hh:122
TpfaFicksLawCache Cache
state the type for the corresponding cache and its filler
Definition: cctpfa/fickslaw.hh:126
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: cctpfa/fickslaw.hh:129
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...