12#ifndef DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
13#define DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
15#include <dune/common/fmatrix.hh>
29template <
class TypeTag,
class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
30class MaxwellStefansLawImplementation;
36template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
42 using FVElementGeometry =
typename GridGeometry::LocalView;
43 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
45 using GridView =
typename GridGeometry::GridView;
46 using Element =
typename GridView::template Codim<0>::Entity;
54 static const int numComponents = ModelTraits::numFluidComponents();
55 static const int numPhases = ModelTraits::numFluidPhases();
56 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
58 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
59 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
61 static_assert(ModelTraits::numFluidPhases() == 1,
"Only one phase allowed supported!");
72 {
return referenceSystem; }
86 template<
class ElementVolumeVariables>
87 static CellCenterPrimaryVariables
flux(
const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace& scvf)
95 CellCenterPrimaryVariables componentFlux(0.0);
96 ReducedComponentVector moleFracInside(0.0);
97 ReducedComponentVector moleFracOutside(0.0);
98 ReducedComponentVector reducedFlux(0.0);
99 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
100 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
103 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
104 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
105 const auto rhoInside = insideVolVars.density();
106 const auto rhoOutside = outsideVolVars.density();
109 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
111 if(compIdx == FluidSystem::getMainComponent(0))
115 const auto eqIdx = Indices::conti0EqIdx + compIdx;
118 const auto bcTypes = problem.boundaryTypes(element, scvf);
119 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
120 return componentFlux;
125 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
128 const auto xInside = insideVolVars.moleFraction(compIdx);
130 const auto xOutside = outsideVolVars.moleFraction(compIdx);
132 moleFracInside[compIdx] = xInside;
133 moleFracOutside[compIdx] = xOutside;
137 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
139 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
141 const auto insideScvIdx = scvf.insideScvIdx();
142 const auto& insideScv = fvGeometry.scv(insideScvIdx);
143 const auto outsideScvIdx = scvf.outsideScvIdx();
144 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
146 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
148 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
153 moleFracOutside -= moleFracInside;
154 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
155 reducedFlux *= omegai*rhoInside;
159 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
160 const Scalar omegaj = calculateOmega_(outsideDistance, outsideVolVars.extrusionFactor());
162 reducedDiffusionMatrixInside.invert();
163 reducedDiffusionMatrixOutside.invert();
164 reducedDiffusionMatrixInside *= omegai*rhoInside;
165 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
168 ReducedComponentVector helperVector(0.0);
169 ReducedComponentVector gradientVectori(0.0);
170 ReducedComponentVector gradientVectorj(0.0);
172 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
173 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
175 auto gradientVectorij = (gradientVectori + gradientVectorj);
178 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
180 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
183 helperVector -=moleFracInside;
184 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
187 reducedFlux *= -Extrusion::area(fvGeometry, scvf);
189 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
191 componentFlux[compIdx] = reducedFlux[compIdx];
192 componentFlux[numComponents-1] -= reducedFlux[compIdx];
195 return componentFlux ;
199 static Scalar calculateOmega_(
const Scalar
distance,
200 const Scalar extrusionFactor)
203 omega *= extrusionFactor;
208 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
209 const FVElementGeometry& fvGeometry,
210 const VolumeVariables& volVars,
211 const SubControlVolumeFace& scvf)
213 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
215 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
217 const auto xi = volVars.moleFraction(compIIdx);
218 const auto avgMolarMass = volVars.averageMolarMass(0);
219 const auto Mn = FluidSystem::molarMass(numComponents-1);
220 const Scalar tin = getEffectiveDiffusionCoefficient_(volVars, compIIdx, numComponents-1);
223 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
225 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
228 if (compJIdx == compIIdx)
231 const auto xj = volVars.moleFraction(compJIdx);
232 const auto Mi = FluidSystem::molarMass(compIIdx);
233 const auto Mj = FluidSystem::molarMass(compJIdx);
234 const Scalar tij = getEffectiveDiffusionCoefficient_(volVars, compIIdx, compJIdx);
235 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
236 if (compJIdx < numComponents-1)
237 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
240 return reducedDiffusionMatrix;
243 static Scalar getEffectiveDiffusionCoefficient_(
const VolumeVariables& volVars,
const int phaseIdx,
const int compIdx)
245 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
Definition: maxwellstefandiffusioncoefficients.hh:33
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/maxwellstefanslaw.hh:71
static CellCenterPrimaryVariables flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: staggered/freeflow/maxwellstefanslaw.hh:87
Definition: box/maxwellstefanslaw.hh:33
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Classes related to flux variables caching.
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
The available discretization methods in Dumux.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: method.hh:114
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:20
Definition: fluxvariablescaching.hh:55