25#ifndef DUMUX_DISCRETIZATION_CC_TPFA_MAXWELL_STEFAN_LAW_HH
26#define DUMUX_DISCRETIZATION_CC_TPFA_MAXWELL_STEFAN_LAW_HH
28#include <dune/common/fmatrix.hh>
29#include <dune/common/float_cmp.hh>
43template <
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
44class MaxwellStefansLawImplementation;
50template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
56 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolume =
typename GridGeometry::SubControlVolume;
58 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
60 using GridView =
typename GridGeometry::GridView;
63 using Element =
typename GridView::template Codim<0>::Entity;
65 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
66 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
69 static const int dim = GridView::dimension;
70 static const int dimWorld = GridView::dimensionworld;
77 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
78 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
79 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
88 {
return referenceSystem; }
102 static ComponentFluxVector
flux(
const Problem& problem,
103 const Element& element,
104 const FVElementGeometry& fvGeometry,
105 const ElementVolumeVariables& elemVolVars,
106 const SubControlVolumeFace& scvf,
108 const ElementFluxVariablesCache& elemFluxVarsCache)
112 ComponentFluxVector componentFlux(0.0);
113 ReducedComponentVector moleFracInside(0.0);
114 ReducedComponentVector moleFracOutside(0.0);
115 ReducedComponentVector reducedFlux(0.0);
116 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
117 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
120 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
121 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
122 const auto rhoInside = insideVolVars.density(phaseIdx);
123 const auto rhoOutside = outsideVolVars.density(phaseIdx);
125 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
128 const auto xInside = insideVolVars.moleFraction(phaseIdx, compIdx);
130 const auto xOutside = outsideVolVars.moleFraction(phaseIdx, compIdx);
132 moleFracInside[compIdx] = xInside;
133 moleFracOutside[compIdx] = xOutside;
137 if(!(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0)))
139 const auto insideScvIdx = scvf.insideScvIdx();
140 const auto& insideScv = fvGeometry.scv(insideScvIdx);
141 const Scalar omegai = calculateOmega_(scvf,
143 insideVolVars.extrusionFactor());
146 reducedDiffusionMatrixInside = setupMSMatrix_(problem, element, fvGeometry, insideVolVars, insideScv, phaseIdx);
149 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
151 moleFracOutside -= moleFracInside;
152 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
153 reducedFlux *= omegai*rhoInside;
158 const auto outsideScvIdx = scvf.outsideScvIdx();
159 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
161 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, element, fvGeometry, outsideVolVars, outsideScv, phaseIdx);
165 omegaj = -1.0*calculateOmega_(scvf,
167 outsideVolVars.extrusionFactor());
169 omegaj = calculateOmega_(fvGeometry.flipScvf(scvf.index()),
171 outsideVolVars.extrusionFactor());
173 reducedDiffusionMatrixInside.invert();
174 reducedDiffusionMatrixOutside.invert();
175 reducedDiffusionMatrixInside *= omegai*rhoInside;
176 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
179 ReducedComponentVector helperVector(0.0);
180 ReducedComponentVector gradientVectori(0.0);
181 ReducedComponentVector gradientVectorj(0.0);
183 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
184 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
186 auto gradientVectorij = (gradientVectori + gradientVectorj);
189 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
191 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
194 helperVector -=moleFracInside;
195 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
198 reducedFlux *= -Extrusion::area(fvGeometry, scvf);
199 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
201 componentFlux[compIdx] = reducedFlux[compIdx];
202 componentFlux[numComponents-1] -=reducedFlux[compIdx];
205 return componentFlux ;
209 static Scalar calculateOmega_(
const SubControlVolumeFace& scvf,
210 const SubControlVolume &scv,
211 const Scalar extrusionFactor)
213 auto distanceVector = scvf.ipGlobal();
214 distanceVector -= scv.center();
215 distanceVector /= distanceVector.two_norm2();
217 Scalar omega = (distanceVector * scvf.unitOuterNormal());
218 omega *= extrusionFactor;
223 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
224 const Element& element,
225 const FVElementGeometry& fvGeometry,
226 const VolumeVariables& volVars,
227 const SubControlVolume& scv,
230 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
233 if(Dune::FloatCmp::eq<Scalar>(volVars.saturation(phaseIdx), 0))
234 return reducedDiffusionMatrix;
236 const auto avgMolarMass = volVars.averageMolarMass(phaseIdx);
238 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
240 const auto xi = volVars.moleFraction(phaseIdx, compIIdx);
242 const auto Mn = FluidSystem::molarMass(numComponents-1);
243 Scalar tin = volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1);
246 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
249 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
252 if (compJIdx == compIIdx)
255 const auto xj = volVars.moleFraction(phaseIdx, compJIdx);
256 const auto Mi = FluidSystem::molarMass(compIIdx);
257 const auto Mj = FluidSystem::molarMass(compJIdx);
258 Scalar tij = volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
259 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
261 if (compJIdx < numComponents-1)
262 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
265 return reducedDiffusionMatrix;
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 the Maxwell- Stefan diffusion law....
Classes related to flux variables caching.
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
Definition: box/maxwellstefanslaw.hh:45
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: cctpfa/maxwellstefanslaw.hh:87
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: cctpfa/maxwellstefanslaw.hh:102
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
Definition: fluxvariablescaching.hh:67
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
Definition: maxwellstefandiffusioncoefficients.hh:45
Declares all properties used in Dumux.