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;
68 static const int dim = GridView::dimension;
69 static const int dimWorld = GridView::dimensionworld;
76 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
77 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
78 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
87 {
return referenceSystem; }
101 static ComponentFluxVector
flux(
const Problem& problem,
102 const Element& element,
103 const FVElementGeometry& fvGeometry,
104 const ElementVolumeVariables& elemVolVars,
105 const SubControlVolumeFace& scvf,
107 const ElementFluxVariablesCache& elemFluxVarsCache)
111 ComponentFluxVector componentFlux(0.0);
112 ReducedComponentVector moleFracInside(0.0);
113 ReducedComponentVector moleFracOutside(0.0);
114 ReducedComponentVector reducedFlux(0.0);
115 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
116 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
119 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
120 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
121 const auto rhoInside = insideVolVars.density(phaseIdx);
122 const auto rhoOutside = outsideVolVars.density(phaseIdx);
124 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
127 const auto xInside = insideVolVars.moleFraction(phaseIdx, compIdx);
129 const auto xOutside = outsideVolVars.moleFraction(phaseIdx, compIdx);
131 moleFracInside[compIdx] = xInside;
132 moleFracOutside[compIdx] = xOutside;
136 if(!(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0)))
138 const auto insideScvIdx = scvf.insideScvIdx();
139 const auto& insideScv = fvGeometry.scv(insideScvIdx);
140 const Scalar omegai = calculateOmega_(scvf,
142 insideVolVars.extrusionFactor());
145 reducedDiffusionMatrixInside = setupMSMatrix_(problem, element, fvGeometry, insideVolVars, insideScv, phaseIdx);
148 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
150 moleFracOutside -= moleFracInside;
151 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
152 reducedFlux *= omegai*rhoInside;
157 const auto outsideScvIdx = scvf.outsideScvIdx();
158 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
160 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, element, fvGeometry, outsideVolVars, outsideScv, phaseIdx);
164 omegaj = -1.0*calculateOmega_(scvf,
166 outsideVolVars.extrusionFactor());
168 omegaj = calculateOmega_(fvGeometry.flipScvf(scvf.index()),
170 outsideVolVars.extrusionFactor());
172 reducedDiffusionMatrixInside.invert();
173 reducedDiffusionMatrixOutside.invert();
174 reducedDiffusionMatrixInside *= omegai*rhoInside;
175 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
178 ReducedComponentVector helperVector(0.0);
179 ReducedComponentVector gradientVectori(0.0);
180 ReducedComponentVector gradientVectorj(0.0);
182 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
183 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
185 auto gradientVectorij = (gradientVectori + gradientVectorj);
188 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
190 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
193 helperVector -=moleFracInside;
194 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
197 reducedFlux *= -Extrusion::area(scvf);
198 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
200 componentFlux[compIdx] = reducedFlux[compIdx];
201 componentFlux[numComponents-1] -=reducedFlux[compIdx];
204 return componentFlux ;
208 static Scalar calculateOmega_(
const SubControlVolumeFace& scvf,
209 const SubControlVolume &scv,
210 const Scalar extrusionFactor)
212 auto distanceVector = scvf.ipGlobal();
213 distanceVector -= scv.center();
214 distanceVector /= distanceVector.two_norm2();
216 Scalar omega = (distanceVector * scvf.unitOuterNormal());
217 omega *= extrusionFactor;
222 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
223 const Element& element,
224 const FVElementGeometry& fvGeometry,
225 const VolumeVariables& volVars,
226 const SubControlVolume& scv,
229 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
232 if(Dune::FloatCmp::eq<Scalar>(volVars.saturation(phaseIdx), 0))
233 return reducedDiffusionMatrix;
235 const auto avgMolarMass = volVars.averageMolarMass(phaseIdx);
237 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
239 const auto xi = volVars.moleFraction(phaseIdx, compIIdx);
241 const auto Mn = FluidSystem::molarMass(numComponents-1);
242 Scalar tin = volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1);
245 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
248 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
251 if (compJIdx == compIIdx)
254 const auto xj = volVars.moleFraction(phaseIdx, compJIdx);
255 const auto Mi = FluidSystem::molarMass(compIIdx);
256 const auto Mj = FluidSystem::molarMass(compJIdx);
257 Scalar tij = volVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
258 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
260 if (compJIdx < numComponents-1)
261 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
264 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
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
Definition: propertysystem.hh:150
Definition: box/maxwellstefanslaw.hh:45
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: cctpfa/maxwellstefanslaw.hh:86
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:101
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.