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,
class DiscretizationMethod, 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!");
84 {
return referenceSystem; }
98 template<
class ElementVolumeVariables>
99 static CellCenterPrimaryVariables
flux(
const Problem& problem,
100 const Element& element,
101 const FVElementGeometry& fvGeometry,
102 const ElementVolumeVariables& elemVolVars,
103 const SubControlVolumeFace& scvf)
107 CellCenterPrimaryVariables componentFlux(0.0);
108 ReducedComponentVector moleFracInside(0.0);
109 ReducedComponentVector moleFracOutside(0.0);
110 ReducedComponentVector reducedFlux(0.0);
111 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
112 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
115 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
116 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
117 const auto rhoInside = insideVolVars.density();
118 const auto rhoOutside = outsideVolVars.density();
121 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
123 if(compIdx == FluidSystem::getMainComponent(0))
127 const auto eqIdx = Indices::conti0EqIdx + compIdx;
130 const auto bcTypes = problem.boundaryTypes(element, scvf);
131 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
132 return componentFlux;
137 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
140 const auto xInside = insideVolVars.moleFraction(compIdx);
142 const auto xOutside = outsideVolVars.moleFraction(compIdx);
144 moleFracInside[compIdx] = xInside;
145 moleFracOutside[compIdx] = xOutside;
149 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
151 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
153 const auto insideScvIdx = scvf.insideScvIdx();
154 const auto& insideScv = fvGeometry.scv(insideScvIdx);
155 const auto outsideScvIdx = scvf.outsideScvIdx();
156 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
158 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
160 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
165 moleFracOutside -= moleFracInside;
166 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
167 reducedFlux *= omegai*rhoInside;
171 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
172 const Scalar omegaj = calculateOmega_(outsideDistance, outsideVolVars.extrusionFactor());
174 reducedDiffusionMatrixInside.invert();
175 reducedDiffusionMatrixOutside.invert();
176 reducedDiffusionMatrixInside *= omegai*rhoInside;
177 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
180 ReducedComponentVector helperVector(0.0);
181 ReducedComponentVector gradientVectori(0.0);
182 ReducedComponentVector gradientVectorj(0.0);
184 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
185 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
187 auto gradientVectorij = (gradientVectori + gradientVectorj);
190 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
192 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
195 helperVector -=moleFracInside;
196 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
199 reducedFlux *= -Extrusion::area(fvGeometry, scvf);
201 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
203 componentFlux[compIdx] = reducedFlux[compIdx];
204 componentFlux[numComponents-1] -= reducedFlux[compIdx];
207 return componentFlux ;
211 static Scalar calculateOmega_(
const Scalar
distance,
212 const Scalar extrusionFactor)
215 omega *= extrusionFactor;
220 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
221 const FVElementGeometry& fvGeometry,
222 const VolumeVariables& volVars,
223 const SubControlVolumeFace& scvf)
225 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
227 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
229 const auto xi = volVars.moleFraction(compIIdx);
230 const auto avgMolarMass = volVars.averageMolarMass(0);
231 const auto Mn = FluidSystem::molarMass(numComponents-1);
232 const Scalar tin = getEffectiveDiffusionCoefficient_(volVars, compIIdx, numComponents-1);
235 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
237 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
240 if (compJIdx == compIIdx)
243 const auto xj = volVars.moleFraction(compJIdx);
244 const auto Mi = FluidSystem::molarMass(compIIdx);
245 const auto Mj = FluidSystem::molarMass(compJIdx);
246 const Scalar tij = getEffectiveDiffusionCoefficient_(volVars, compIIdx, compJIdx);
247 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
248 if (compJIdx < numComponents-1)
249 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
252 return reducedDiffusionMatrix;
255 static Scalar getEffectiveDiffusionCoefficient_(
const VolumeVariables& volVars,
const int phaseIdx,
const int compIdx)
257 return volVars.effectiveDiffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
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.
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:294
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Definition: method.hh:103
Definition: box/maxwellstefanslaw.hh:45
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
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/maxwellstefanslaw.hh:83
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:99
Declares all properties used in Dumux.