25#ifndef DUMUX_DISCRETIZATION_CC_TPFA_MAXWELL_STEFAN_LAW_HH
26#define DUMUX_DISCRETIZATION_CC_TPFA_MAXWELL_STEFAN_LAW_HH
28#include <dune/common/float_cmp.hh>
40template <
class TypeTag, DiscretizationMethod discMethod, ReferenceSystemFormulation referenceSystem>
41class MaxwellStefansLawImplementation;
47template <
class TypeTag, ReferenceSystemFormulation referenceSystem>
54 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
55 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
59 using Element =
typename GridView::template Codim<0>::Entity;
64 static const int dim = GridView::dimension;
65 static const int dimWorld = GridView::dimensionworld;
72 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
73 using ReducedComponentVector = Dune::FieldVector<Scalar, numComponents-1>;
74 using ReducedComponentMatrix = Dune::FieldMatrix<Scalar, numComponents-1, numComponents-1>;
81 {
return referenceSystem; }
88 static ComponentFluxVector
flux(
const Problem& problem,
89 const Element& element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const SubControlVolumeFace& scvf,
94 const ElementFluxVariablesCache& elemFluxVarsCache)
98 ComponentFluxVector componentFlux(0.0);
99 ReducedComponentVector moleFracInside(0.0);
100 ReducedComponentVector moleFracOutside(0.0);
101 ReducedComponentVector reducedFlux(0.0);
102 ReducedComponentMatrix reducedDiffusionMatrixInside(0.0);
103 ReducedComponentMatrix reducedDiffusionMatrixOutside(0.0);
106 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
107 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
108 const auto rhoInside = insideVolVars.density(phaseIdx);
109 const auto rhoOutside = outsideVolVars.density(phaseIdx);
111 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
114 const auto xInside = insideVolVars.moleFraction(phaseIdx, compIdx);
116 const auto xOutside = outsideVolVars.moleFraction(phaseIdx, compIdx);
118 moleFracInside[compIdx] = xInside;
119 moleFracOutside[compIdx] = xOutside;
123 if(!(Dune::FloatCmp::eq<Scalar>(insideVolVars.saturation(phaseIdx), 0) || Dune::FloatCmp::eq<Scalar>(outsideVolVars.saturation(phaseIdx), 0)))
125 const auto insideScvIdx = scvf.insideScvIdx();
126 const auto& insideScv = fvGeometry.scv(insideScvIdx);
127 const Scalar omegai = calculateOmega_(scvf,
129 insideVolVars.extrusionFactor());
132 reducedDiffusionMatrixInside = setupMSMatrix_(problem, element, fvGeometry, insideVolVars, insideScv, phaseIdx);
135 if (scvf.boundary() || scvf.numOutsideScvs() > 1)
137 moleFracOutside -= moleFracInside;
138 reducedDiffusionMatrixInside.solve(reducedFlux, moleFracOutside);
139 reducedFlux *= omegai*rhoInside;
144 const auto outsideScvIdx = scvf.outsideScvIdx();
145 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
147 reducedDiffusionMatrixOutside = setupMSMatrix_(problem, element, fvGeometry, outsideVolVars, outsideScv, phaseIdx);
151 omegaj = -1.0*calculateOmega_(scvf,
153 outsideVolVars.extrusionFactor());
155 omegaj = calculateOmega_(fvGeometry.flipScvf(scvf.index()),
157 outsideVolVars.extrusionFactor());
159 reducedDiffusionMatrixInside.invert();
160 reducedDiffusionMatrixOutside.invert();
161 reducedDiffusionMatrixInside *= omegai*rhoInside;
162 reducedDiffusionMatrixOutside *= omegaj*rhoOutside;
165 ReducedComponentVector helperVector(0.0);
166 ReducedComponentVector gradientVectori(0.0);
167 ReducedComponentVector gradientVectorj(0.0);
169 reducedDiffusionMatrixInside.mv(moleFracInside, gradientVectori);
170 reducedDiffusionMatrixOutside.mv(moleFracOutside, gradientVectorj);
172 auto gradientVectorij = (gradientVectori + gradientVectorj);
175 reducedDiffusionMatrixOutside += reducedDiffusionMatrixInside;
177 reducedDiffusionMatrixOutside.solve(helperVector, gradientVectorij);
180 helperVector -=moleFracInside;
181 reducedDiffusionMatrixInside.mv(helperVector, reducedFlux);
184 reducedFlux *= -scvf.area();
185 for (
int compIdx = 0; compIdx < numComponents-1; compIdx++)
187 componentFlux[compIdx] = reducedFlux[compIdx];
188 componentFlux[numComponents-1] -=reducedFlux[compIdx];
191 return componentFlux ;
195 static Scalar calculateOmega_(
const SubControlVolumeFace& scvf,
196 const SubControlVolume &scv,
197 const Scalar extrusionFactor)
199 auto distanceVector = scvf.ipGlobal();
200 distanceVector -= scv.center();
201 distanceVector /= distanceVector.two_norm2();
203 Scalar omega = (distanceVector * scvf.unitOuterNormal());
204 omega *= extrusionFactor;
209 static ReducedComponentMatrix setupMSMatrix_(
const Problem& problem,
210 const Element& element,
211 const FVElementGeometry& fvGeometry,
212 const VolumeVariables& volVars,
213 const SubControlVolume& scv,
216 ReducedComponentMatrix reducedDiffusionMatrix(0.0);
219 if(Dune::FloatCmp::eq<Scalar>(volVars.saturation(phaseIdx), 0))
220 return reducedDiffusionMatrix;
222 const auto avgMolarMass = volVars.averageMolarMass(phaseIdx);
224 for (
int compIIdx = 0; compIIdx < numComponents-1; compIIdx++)
226 const auto xi = volVars.moleFraction(phaseIdx, compIIdx);
228 const auto Mn = FluidSystem::molarMass(numComponents-1);
229 Scalar tin = getDiffusionCoefficient(phaseIdx, compIIdx, numComponents-1, problem, element, volVars, scv);
230 tin = EffDiffModel::effectiveDiffusivity(volVars.porosity(), volVars.saturation(phaseIdx), tin);
233 reducedDiffusionMatrix[compIIdx][compIIdx] += xi*avgMolarMass/(tin*Mn);
236 for (
int compJIdx = 0; compJIdx < numComponents; compJIdx++)
239 if (compJIdx == compIIdx)
242 const auto xj = volVars.moleFraction(phaseIdx, compJIdx);
243 const auto Mi = FluidSystem::molarMass(compIIdx);
244 const auto Mj = FluidSystem::molarMass(compJIdx);
245 Scalar tij = getDiffusionCoefficient(phaseIdx, compIIdx, compJIdx, problem, element, volVars, scv);
246 tij = EffDiffModel::effectiveDiffusivity(volVars.porosity(), volVars.saturation(phaseIdx), tij);
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 template <class T = TypeTag, typename std::enable_if_t<GetPropType<T, Properties::FluidSystem>::isTracerFluidSystem(),
int> =0 >
256 static Scalar getDiffusionCoefficient(
const int phaseIdx,
259 const Problem& problem,
260 const Element& element,
261 const VolumeVariables& volVars,
262 const SubControlVolume& scv)
264 return FluidSystem::binaryDiffusionCoefficient(compIIdx,
271 template <class T = TypeTag, typename std::enable_if_t<!GetPropType<T, Properties::FluidSystem>::isTracerFluidSystem(),
int> =0 >
272 static Scalar getDiffusionCoefficient(
const int phaseIdx,
275 const Problem& problem,
276 const Element& element,
277 const VolumeVariables& volVars,
278 const SubControlVolume& scv)
280 auto fluidState = volVars.fluidState();
281 typename FluidSystem::ParameterCache paramCache;
282 paramCache.updateAll(fluidState);
283 return FluidSystem::binaryDiffusionCoefficient(fluidState,
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
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: cctpfa/maxwellstefanslaw.hh:80
static ComponentFluxVector flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Definition: cctpfa/maxwellstefanslaw.hh:88
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
Definition: fluxvariablescaching.hh:64
Declares all properties used in Dumux.