26#ifndef DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
27#define DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
37template<
class TypeTag,
bool enableChemicalNonEquilibrium>
40template <
class TypeTag>
48template<
class TypeTag>
55 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
57 using Element =
typename GridView::template Codim<0>::Entity;
63 using ParentType::ParentType;
75 const Element& element,
76 const FVElementGeometry& fvGeometry,
77 const ElementVolumeVariables& elemVolVars,
78 const SubControlVolume &scv)
const
80 NumEqVector source(0.0);
83 source += problem.source(element, fvGeometry, elemVolVars, scv);
86 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
90 if constexpr(ModelTraits::enableThermalNonEquilibrium())
92 EnergyLocalResidual::computeSourceEnergy(source,
108template<
class TypeTag>
115 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
116 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
121 using Element =
typename GridView::template Codim<0>::Entity;
128 using Indices =
typename ModelTraits::Indices;
130 static constexpr int numPhases = ModelTraits::numFluidPhases();
131 static constexpr int numComponents = ModelTraits::numFluidComponents();
133 static constexpr auto conti0EqIdx = Indices::conti0EqIdx;
134 static constexpr auto comp0Idx = FluidSystem::comp0Idx;
135 static constexpr auto phase0Idx = FluidSystem::phase0Idx;
136 static constexpr auto phase1Idx = FluidSystem::phase1Idx;
138 static_assert(numPhases > 1,
139 "chemical non-equlibrium only makes sense for multiple phases");
140 static_assert(numPhases == numComponents,
141 "currently chemical non-equilibrium is only available when numPhases equals numComponents");
142 static_assert(ModelTraits::useMoles(),
143 "chemical nonequilibrium can only be calculated based on mole fractions not mass fractions");
146 using ParentType::ParentType;
155 const SubControlVolume& scv,
156 const VolumeVariables& volVars)
const
158 NumEqVector storage(0.0);
161 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
163 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
165 auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
166 storage[eqIdx] += volVars.porosity()
167 * volVars.saturation(phaseIdx)
168 * volVars.molarDensity(phaseIdx)
169 * volVars.moleFraction(phaseIdx, compIdx);
172 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
175 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
190 const Element& element,
191 const FVElementGeometry& fvGeometry,
192 const ElementVolumeVariables& elemVolVars,
193 const SubControlVolumeFace& scvf,
194 const ElementFluxVariablesCache& elemFluxVarsCache)
const
196 FluxVariables fluxVars;
197 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
198 NumEqVector flux(0.0);
201 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
203 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
204 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
207 const auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
210 const auto upwindTerm = [phaseIdx, compIdx] (
const auto& volVars)
211 {
return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
213 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
216 if (compIdx == phaseIdx)
220 flux[eqIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
222 flux[eqIdx] += diffusiveFluxes[compIdx];
225 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
229 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
244 const Element& element,
245 const FVElementGeometry& fvGeometry,
246 const ElementVolumeVariables& elemVolVars,
247 const SubControlVolume &scv)
const
249 NumEqVector source(0.0);
253 const auto& volVars = elemVolVars[scv];
254 std::array<std::array<Scalar, numComponents>, numPhases> componentIntoPhaseMassTransfer = {{{0.0},{0.0}}};
257 const Scalar characteristicLength = volVars.characteristicLength() ;
258 const Scalar factorMassTransfer = volVars.factorMassTransfer() ;
260 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
262 for(
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
264 const Scalar sherwoodNumber = volVars.sherwoodNumber(phaseIdx);
266 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
268 if (compIdx <= numPhases)
271 if (phaseIdx == compIdx)
274 const Scalar xNonEquil = volVars.moleFraction(phaseIdx, compIdx);
277 const Scalar xEquil = volVars.xEquil(phaseIdx, compIdx);
279 const Scalar diffCoeff = volVars.diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
282 const Scalar compFluxIntoOtherPhase = factorMassTransfer * (xEquil-xNonEquil)/characteristicLength * awn * volVars.molarDensity(phaseIdx) * diffCoeff * sherwoodNumber;
284 componentIntoPhaseMassTransfer[phaseIdx][compIdx] += compFluxIntoOtherPhase;
285 componentIntoPhaseMassTransfer[compIdx][compIdx] -= compFluxIntoOtherPhase;
292 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
294 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
296 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
297 source[eqIdx] += componentIntoPhaseMassTransfer[phaseIdx][compIdx];
300 if (!isfinite(source[eqIdx]))
305 if constexpr (ModelTraits::enableThermalNonEquilibrium())
309 EnergyLocalResidual::computeSourceEnergy(source,
317 source += problem.source(element, fvGeometry, elemVolVars, scv);
320 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
A helper to deduce a vector with the same size as numbers of equations.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Definition: porousmediumflow/nonequilibrium/localresidual.hh:38
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculates the source term of the equation.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:74
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculates the source term of the equation.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:243
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Calculates the storage for all mass balance equations.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:154
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Calculates the flux for all mass balance equations.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:189
Declares all properties used in Dumux.
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...