24#ifndef DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
27#include <dune/common/fmatrix.hh>
41template <
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
42class MaxwellStefansLawImplementation;
48template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
54 using FVElementGeometry =
typename GridGeometry::LocalView;
55 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
56 using GridView =
typename GridGeometry::GridView;
57 using Element =
typename GridView::template Codim<0>::Entity;
65 static const int numComponents = ModelTraits::numFluidComponents();
66 static const int numPhases = ModelTraits::numFluidPhases();
67 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
69 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
70 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
72 static_assert(ModelTraits::numFluidPhases() == 1,
"Only one phase allowed supported!");
81 {
return referenceSystem; }
90 template<
class ElementVolumeVariables>
91 static CellCenterPrimaryVariables
flux(
const Problem& problem,
92 const Element& element,
93 const FVElementGeometry& fvGeometry,
94 const ElementVolumeVariables& elemVolVars,
95 const SubControlVolumeFace& scvf)
99 CellCenterPrimaryVariables componentFlux(0.0);
100 ReducedComponentVector moleFracInside(0.0);
101 ReducedComponentVector moleFracOutside(0.0);
102 ReducedComponentVector reducedFlux(0.0);
103 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
104 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
107 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
108 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
109 const auto rhoInside = insideVolVars.density();
110 const auto rhoOutside = outsideVolVars.density();
113 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
115 if(compIdx == FluidSystem::getMainComponent(0))
119 const auto eqIdx = Indices::conti0EqIdx + compIdx;
122 const auto bcTypes = problem.boundaryTypes(element, scvf);
123 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
124 return componentFlux;
129 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
132 const auto xInside = insideVolVars.moleFraction(compIdx);
134 const auto xOutside = outsideVolVars.moleFraction(compIdx);
136 moleFracInside[compIdx] = xInside;
137 moleFracOutside[compIdx] = xOutside;
141 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
143 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
145 const auto insideScvIdx = scvf.insideScvIdx();
146 const auto& insideScv = fvGeometry.scv(insideScvIdx);
147 const auto outsideScvIdx = scvf.outsideScvIdx();
148 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
150 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
152 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
157 moleFracOutside -= moleFracInside;
158 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
159 reducedFlux *= omegai*rhoInside;
163 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
164 const Scalar omegaj = calculateOmega_(outsideDistance, outsideVolVars.extrusionFactor());
166 reducedDiffusionMatrixInside.invert();
167 reducedDiffusionMatrixOutside.invert();
168 reducedDiffusionMatrixInside *= omegai*rhoInside;
169 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
172 ReducedComponentVector helperVector(0.0);
173 ReducedComponentVector gradientVectori(0.0);
174 ReducedComponentVector gradientVectorj(0.0);
176 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
177 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
179 auto gradientVectorij = (gradientVectori + gradientVectorj);
182 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
184 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
187 helperVector -=moleFracInside;
188 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
191 reducedFlux *= -scvf.area();
193 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
195 componentFlux[compIdx] = reducedFlux[compIdx];
196 componentFlux[numComponents-1] -= reducedFlux[compIdx];
199 return componentFlux ;
203 static Scalar calculateOmega_(
const Scalar distance,
204 const Scalar extrusionFactor)
206 Scalar omega = 1/distance;
207 omega *= extrusionFactor;
212 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
213 const FVElementGeometry& fvGeometry,
214 const VolumeVariables& volVars,
215 const SubControlVolumeFace& scvf)
217 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
219 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
221 const auto xi = volVars.moleFraction(compIIdx);
222 const auto avgMolarMass = volVars.averageMolarMass(0);
223 const auto Mn = FluidSystem::molarMass(numComponents-1);
224 const Scalar tin = getEffectiveDiffusionCoefficient_(volVars, compIIdx, numComponents-1);
227 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
229 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
232 if (compJIdx == compIIdx)
235 const auto xj = volVars.moleFraction(compJIdx);
236 const auto Mi = FluidSystem::molarMass(compIIdx);
237 const auto Mj = FluidSystem::molarMass(compJIdx);
238 const Scalar tij = getEffectiveDiffusionCoefficient_(volVars, compIIdx, compJIdx);
239 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
240 if (compJIdx < numComponents-1)
241 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
244 return reducedDiffusionMatrix;
247 static Scalar getEffectiveDiffusionCoefficient_(
const VolumeVariables& volVars,
const int phaseIdx,
const int compIdx)
249 if constexpr (Dumux::Deprecated::hasEffDiffCoeff<VolumeVariables>)
250 return volVars.effectiveDiffusionCoefficient(phaseIdx,
251 VolumeVariables::FluidSystem::getMainComponent(phaseIdx), compIdx);
255 return volVars.effectiveDiffusivity(phaseIdx, compIdx);
Classes related to flux variables caching.
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
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
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
static CellCenterPrimaryVariables flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: staggered/freeflow/maxwellstefanslaw.hh:91
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/maxwellstefanslaw.hh:80
Declares all properties used in Dumux.