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 GridGeometry::SubControlVolumeFace;
57 using GridView =
typename GridGeometry::GridView;
58 using Element =
typename GridView::template Codim<0>::Entity;
66 static const int numComponents = ModelTraits::numFluidComponents();
67 static const int numPhases = ModelTraits::numFluidPhases();
68 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
70 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
71 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
73 static_assert(ModelTraits::numFluidPhases() == 1,
"Only one phase allowed supported!");
82 {
return referenceSystem; }
96 template<
class ElementVolumeVariables>
97 static CellCenterPrimaryVariables
flux(
const Problem& problem,
98 const Element& element,
99 const FVElementGeometry& fvGeometry,
100 const ElementVolumeVariables& elemVolVars,
101 const SubControlVolumeFace& scvf)
105 CellCenterPrimaryVariables componentFlux(0.0);
106 ReducedComponentVector moleFracInside(0.0);
107 ReducedComponentVector moleFracOutside(0.0);
108 ReducedComponentVector reducedFlux(0.0);
109 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
110 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
113 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
114 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
115 const auto rhoInside = insideVolVars.density();
116 const auto rhoOutside = outsideVolVars.density();
119 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
121 if(compIdx == FluidSystem::getMainComponent(0))
125 const auto eqIdx = Indices::conti0EqIdx + compIdx;
128 const auto bcTypes = problem.boundaryTypes(element, scvf);
129 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
130 return componentFlux;
135 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
138 const auto xInside = insideVolVars.moleFraction(compIdx);
140 const auto xOutside = outsideVolVars.moleFraction(compIdx);
142 moleFracInside[compIdx] = xInside;
143 moleFracOutside[compIdx] = xOutside;
147 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
149 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
151 const auto insideScvIdx = scvf.insideScvIdx();
152 const auto& insideScv = fvGeometry.scv(insideScvIdx);
153 const auto outsideScvIdx = scvf.outsideScvIdx();
154 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
156 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
158 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
163 moleFracOutside -= moleFracInside;
164 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
165 reducedFlux *= omegai*rhoInside;
169 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
170 const Scalar omegaj = calculateOmega_(outsideDistance, 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);
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 Scalar
distance,
210 const Scalar extrusionFactor)
213 omega *= extrusionFactor;
218 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
219 const FVElementGeometry& fvGeometry,
220 const VolumeVariables& volVars,
221 const SubControlVolumeFace& scvf)
223 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
225 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
227 const auto xi = volVars.moleFraction(compIIdx);
228 const auto avgMolarMass = volVars.averageMolarMass(0);
229 const auto Mn = FluidSystem::molarMass(numComponents-1);
230 const Scalar tin = getEffectiveDiffusionCoefficient_(volVars, compIIdx, numComponents-1);
233 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
235 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
238 if (compJIdx == compIIdx)
241 const auto xj = volVars.moleFraction(compJIdx);
242 const auto Mi = FluidSystem::molarMass(compIIdx);
243 const auto Mj = FluidSystem::molarMass(compJIdx);
244 const Scalar tij = getEffectiveDiffusionCoefficient_(volVars, compIIdx, compJIdx);
245 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
246 if (compJIdx < numComponents-1)
247 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
250 return reducedDiffusionMatrix;
253 static Scalar getEffectiveDiffusionCoefficient_(
const VolumeVariables& volVars,
const int phaseIdx,
const int compIdx)
255 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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....
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: geometry/distance.hh:138
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 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 (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)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: staggered/freeflow/maxwellstefanslaw.hh:97
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/maxwellstefanslaw.hh:81
Declares all properties used in Dumux.