25#ifndef DUMUX_DISCRETIZATION_BOX_MAXWELL_STEFAN_LAW_HH
26#define DUMUX_DISCRETIZATION_BOX_MAXWELL_STEFAN_LAW_HH
28#include <dune/common/float_cmp.hh>
29#include <dune/common/fmatrix.hh>
43template <
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
44class MaxwellStefansLawImplementation;
50template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
58 using FVElementGeometry =
typename GridGeometry::LocalView;
59 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
64 using Element =
typename GridView::template Codim<0>::Entity;
69 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
70 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
71 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
80 {
return referenceSystem; }
87 static ComponentFluxVector
flux(
const Problem& problem,
88 const Element& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const SubControlVolumeFace& scvf,
93 const ElementFluxVariablesCache& elemFluxVarsCache)
97 ComponentFluxVector componentFlux(0.0);
98 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
99 ReducedComponentVector reducedFlux(0.0);
100 ComponentFluxVector moleFrac(0.0);
101 ReducedComponentVector normalX(0.0);
104 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
105 const auto& shapeValues = fluxVarsCache.shapeValues();
108 Scalar avgMolarMass(0.0);
109 for (
auto&& scv : scvs(fvGeometry))
111 const auto& volVars = elemVolVars[scv];
114 rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
116 avgMolarMass += volVars.averageMolarMass(phaseIdx)*shapeValues[scv.indexInElement()][0];
118 for (
int compIdx = 0; compIdx < numComponents; compIdx++)
120 moleFrac[compIdx] += volVars.moleFraction(phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0];
124 reducedDiffusionMatrix = setupMSMatrix_(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, moleFrac, avgMolarMass);
126 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
128 Dune::FieldVector<Scalar, GridView::dimensionworld> gradX(0.0);
129 for (
auto&& scv : scvs(fvGeometry))
131 const auto& volVars = elemVolVars[scv];
134 gradX.axpy(volVars.moleFraction(phaseIdx, compIdx), fluxVarsCache.gradN(scv.indexInElement()));
137 normalX[compIdx] = gradX *scvf.unitOuterNormal();
139 reducedDiffusionMatrix.solve(reducedFlux,normalX);
140 reducedFlux *= -1.0*rho*Extrusion::area(scvf);
142 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
144 componentFlux[compIdx] = reducedFlux[compIdx];
145 componentFlux[numComponents-1] -= reducedFlux[compIdx];
147 return componentFlux ;
152 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
153 const Element& element,
154 const FVElementGeometry& fvGeometry,
155 const ElementVolumeVariables& elemVolVars,
156 const SubControlVolumeFace& scvf,
158 const ComponentFluxVector moleFrac,
159 const Scalar avgMolarMass)
161 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
163 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
164 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
167 if(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0))
168 return reducedDiffusionMatrix;
170 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
173 const auto xi = moleFrac[compIIdx];
174 const auto Mn = FluidSystem::molarMass(numComponents-1);
176 auto tinInside = insideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1);
177 auto tinOutside = outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1);
180 tinInside *= insideVolVars.extrusionFactor();
181 tinOutside *= outsideVolVars.extrusionFactor();
184 const auto tin = problem.spatialParams().harmonicMean(tinInside, tinOutside, scvf.unitOuterNormal());
187 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
190 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
193 if (compIIdx == compJIdx)
197 const auto xj = moleFrac[compJIdx];
198 const auto Mi = FluidSystem::molarMass(compIIdx);
199 const auto Mj = FluidSystem::molarMass(compJIdx);
202 auto tijInside = insideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
203 auto tijOutside = outsideVolVars.effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx);
206 tijInside *= insideVolVars.extrusionFactor();
207 tijOutside *= outsideVolVars.extrusionFactor();
210 const auto tij = problem.spatialParams().harmonicMean(tijInside, tijOutside, scvf.unitOuterNormal());
212 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
213 if (compJIdx < numComponents-1)
214 reducedDiffusionMatrix[compIIdx][compJIdx] +=xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
217 return reducedDiffusionMatrix;
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....
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
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/maxwellstefanslaw.hh:79
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the diffusive fluxes of all components within a fluid phase across the given sub-control volu...
Definition: box/maxwellstefanslaw.hh:87
Container storing the diffusion coefficients required by the Maxwell- Stefan diffusion law....
Definition: maxwellstefandiffusioncoefficients.hh:45
Declares all properties used in Dumux.