24#ifndef DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_MAXWELL_STEFAN_LAW_HH
27#include <dune/common/float_cmp.hh>
39template <
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
40class MaxwellStefansLawImplementation;
46template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
54 using GridView =
typename GridGeometry::GridView;
55 using Element =
typename GridView::template Codim<0>::Entity;
64 static const int numComponents = ModelTraits::numFluidComponents();
65 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
67 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
68 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
70 static_assert(ModelTraits::numFluidPhases() == 1,
"Only one phase allowed supported!");
79 {
return referenceSystem; }
86 static CellCenterPrimaryVariables
flux(
const Problem& problem,
87 const Element& element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const SubControlVolumeFace& scvf)
94 CellCenterPrimaryVariables componentFlux(0.0);
95 ReducedComponentVector moleFracInside(0.0);
96 ReducedComponentVector moleFracOutside(0.0);
97 ReducedComponentVector reducedFlux(0.0);
98 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
99 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
102 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
103 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
104 const auto rhoInside = insideVolVars.density();
105 const auto rhoOutside = outsideVolVars.density();
108 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
110 if(compIdx == FluidSystem::getMainComponent(0))
114 const auto eqIdx = Indices::conti0EqIdx + compIdx;
117 const auto bcTypes = problem.boundaryTypes(element, scvf);
118 if((bcTypes.isOutflow(eqIdx)) || (bcTypes.isSymmetry()))
119 return componentFlux;
124 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
127 const auto xInside = insideVolVars.moleFraction(compIdx);
129 const auto xOutside = outsideVolVars.moleFraction(compIdx);
131 moleFracInside[compIdx] = xInside;
132 moleFracOutside[compIdx] = xOutside;
136 reducedDiffusionMatrixInside = setupMSMatrix_(problem, fvGeometry, insideVolVars, scvf);
138 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, fvGeometry, outsideVolVars, scvf);
140 const auto insideScvIdx = scvf.insideScvIdx();
141 const auto& insideScv = fvGeometry.scv(insideScvIdx);
142 const auto outsideScvIdx = scvf.outsideScvIdx();
143 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
145 const Scalar insideDistance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
147 const Scalar omegai = calculateOmega_(insideDistance, insideVolVars.extrusionFactor());
152 moleFracOutside -= moleFracInside;
153 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
154 reducedFlux *= omegai*rhoInside;
158 const Scalar outsideDistance = (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm();
159 const Scalar omegaj = calculateOmega_(outsideDistance, outsideVolVars.extrusionFactor());
161 reducedDiffusionMatrixInside.invert();
162 reducedDiffusionMatrixOutside.invert();
163 reducedDiffusionMatrixInside *= omegai*rhoInside;
164 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
167 ReducedComponentVector helperVector(0.0);
168 ReducedComponentVector gradientVectori(0.0);
169 ReducedComponentVector gradientVectorj(0.0);
171 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
172 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
174 auto gradientVectorij = (gradientVectori + gradientVectorj);
177 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
179 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
182 helperVector -=moleFracInside;
183 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
186 reducedFlux *= -scvf.area();
188 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
190 componentFlux[compIdx] = reducedFlux[compIdx];
191 componentFlux[numComponents-1] -= reducedFlux[compIdx];
194 return componentFlux ;
198 static Scalar calculateOmega_(
const Scalar distance,
199 const Scalar extrusionFactor)
201 Scalar omega = 1/distance;
202 omega *= extrusionFactor;
207 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
208 const FVElementGeometry& fvGeometry,
209 const VolumeVariables& volVars,
210 const SubControlVolumeFace& scvf)
212 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
214 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
216 const auto xi = volVars.moleFraction(compIIdx);
217 const auto avgMolarMass = volVars.averageMolarMass(0);
218 const auto Mn = FluidSystem::molarMass(numComponents-1);
219 const Scalar tin = volVars.effectiveDiffusivity(compIIdx, numComponents-1);
222 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
224 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
227 if (compJIdx == compIIdx)
230 const auto xj = volVars.moleFraction(compJIdx);
231 const auto Mi = FluidSystem::molarMass(compIIdx);
232 const auto Mj = FluidSystem::molarMass(compJIdx);
233 const Scalar tij = volVars.effectiveDiffusivity(compIIdx, compJIdx);
234 reducedDiffusionMatrix[compIIdx][compIIdx] += xj*avgMolarMass/(tij*Mi);
235 if (compJIdx < numComponents-1)
236 reducedDiffusionMatrix[compIIdx][compJIdx] += xi*(avgMolarMass/(tin*Mn) - avgMolarMass/(tij*Mj));
239 return reducedDiffusionMatrix;
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Classes related to flux variables caching.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: staggered/freeflow/maxwellstefanslaw.hh:78
static CellCenterPrimaryVariables flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: staggered/freeflow/maxwellstefanslaw.hh:86
Declares all properties used in Dumux.