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, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
44class MaxwellStefansLawImplementation;
50template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
57 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
58 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
62 using Element =
typename GridView::template Codim<0>::Entity;
67 static const int dim = GridView::dimension;
68 static const int dimWorld = GridView::dimensionworld;
75 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
76 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
77 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
84 {
return referenceSystem; }
93 static ComponentFluxVector
flux(
const Problem& problem,
94 const Element& element,
95 const FVElementGeometry& fvGeometry,
96 const ElementVolumeVariables& elemVolVars,
97 const SubControlVolumeFace& scvf,
99 const ElementFluxVariablesCache& elemFluxVarsCache)
103 ComponentFluxVector componentFlux(0.0);
104 ReducedComponentVector moleFracInside(0.0);
105 ReducedComponentVector moleFracOutside(0.0);
106 ReducedComponentVector reducedFlux(0.0);
107 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
108 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
111 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
112 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
113 const auto rhoInside = insideVolVars.density(phaseIdx);
114 const auto rhoOutside = outsideVolVars.density(phaseIdx);
116 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
119 const auto xInside = insideVolVars.moleFraction(phaseIdx, compIdx);
121 const auto xOutside = outsideVolVars.moleFraction(phaseIdx, compIdx);
123 moleFracInside[compIdx] = xInside;
124 moleFracOutside[compIdx] = xOutside;
128 if(!(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0)))
130 const auto insideScvIdx = scvf.insideScvIdx();
131 const auto& insideScv = fvGeometry.scv(insideScvIdx);
132 const Scalar omegai = calculateOmega_(scvf,
134 insideVolVars.extrusionFactor());
137 reducedDiffusionMatrixInside = setupMSMatrix_(problem, element, fvGeometry, insideVolVars, insideScv, phaseIdx);
140 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
142 moleFracOutside -= moleFracInside;
143 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
144 reducedFlux *= omegai*rhoInside;
149 const auto outsideScvIdx = scvf.outsideScvIdx();
150 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
152 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, element, fvGeometry, outsideVolVars, outsideScv, phaseIdx);
156 omegaj = -1.0*calculateOmega_(scvf,
158 outsideVolVars.extrusionFactor());
160 omegaj = calculateOmega_(fvGeometry.flipScvf(scvf.index()),
162 outsideVolVars.extrusionFactor());
164 reducedDiffusionMatrixInside.invert();
165 reducedDiffusionMatrixOutside.invert();
166 reducedDiffusionMatrixInside *= omegai*rhoInside;
167 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
170 ReducedComponentVector helperVector(0.0);
171 ReducedComponentVector gradientVectori(0.0);
172 ReducedComponentVector gradientVectorj(0.0);
174 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
175 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
177 auto gradientVectorij = (gradientVectori + gradientVectorj);
180 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
182 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
185 helperVector -=moleFracInside;
186 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
189 reducedFlux *= -scvf.area();
190 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
192 componentFlux[compIdx] = reducedFlux[compIdx];
193 componentFlux[numComponents-1] -=reducedFlux[compIdx];
196 return componentFlux ;
200 static Scalar calculateOmega_(
const SubControlVolumeFace& scvf,
201 const SubControlVolume &scv,
202 const Scalar extrusionFactor)
204 auto distanceVector = scvf.ipGlobal();
205 distanceVector -= scv.center();
206 distanceVector /= distanceVector.two_norm2();
208 Scalar omega = (distanceVector * scvf.unitOuterNormal());
209 omega *= extrusionFactor;
214 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
215 const Element& element,
216 const FVElementGeometry& fvGeometry,
217 const VolumeVariables& volVars,
218 const SubControlVolume& scv,
221 using EffDiffModel = GetPropType<TypeTag, Properties::EffectiveDiffusivityModel>;
222 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
225 if(Dune::FloatCmp::eq<Scalar>(volVars.saturation(phaseIdx), 0))
226 return reducedDiffusionMatrix;
228 const auto avgMolarMass = volVars.averageMolarMass(phaseIdx);
230 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
232 const auto xi = volVars.moleFraction(phaseIdx, compIIdx);
234 const auto Mn = FluidSystem::molarMass(numComponents-1);
235 Scalar tin = Deprecated::template effectiveMSDiffusionCoefficient<EffDiffModel, FluidSystem>(volVars, phaseIdx, compIIdx, numComponents-1, problem, element, scv);
238 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
241 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
244 if (compJIdx == compIIdx)
247 const auto xj = volVars.moleFraction(phaseIdx, compJIdx);
248 const auto Mi = FluidSystem::molarMass(compIIdx);
249 const auto Mj = FluidSystem::molarMass(compJIdx);
250 Scalar tij = Deprecated::template effectiveMSDiffusionCoefficient<EffDiffModel, FluidSystem>(volVars, phaseIdx, compIIdx, compJIdx, problem, element, scv);
251 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
253 if (compJIdx < numComponents-1)
254 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
257 return reducedDiffusionMatrix;
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...
Classes related to flux variables caching.
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
Definition: maxwellstefanslaw.hh:37
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: cctpfa/maxwellstefanslaw.hh:83
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: cctpfa/maxwellstefanslaw.hh:93
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
Definition: fluxvariablescaching.hh:64
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
Definition: maxwellstefandiffusioncoefficients.hh:45
Declares all properties used in Dumux.