26#ifndef DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
27#define DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
36template<
class TypeTag,
bool enableChemicalNonEquilibrium>
39template <
class TypeTag>
47template<
class TypeTag>
54 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
56 using Element =
typename GridView::template Codim<0>::Entity;
62 using ParentType::ParentType;
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const SubControlVolume &scv)
const
79 NumEqVector source(0.0);
82 source += problem.source(element, fvGeometry, elemVolVars, scv);
85 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
89 if constexpr(ModelTraits::enableThermalNonEquilibrium())
91 EnergyLocalResidual::computeSourceEnergy(source,
107template<
class TypeTag>
114 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
115 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
120 using Element =
typename GridView::template Codim<0>::Entity;
127 using Indices =
typename ModelTraits::Indices;
129 static constexpr int numPhases = ModelTraits::numFluidPhases();
130 static constexpr int numComponents = ModelTraits::numFluidComponents();
132 static constexpr auto conti0EqIdx = Indices::conti0EqIdx;
133 static constexpr auto comp0Idx = FluidSystem::comp0Idx;
134 static constexpr auto phase0Idx = FluidSystem::phase0Idx;
135 static constexpr auto phase1Idx = FluidSystem::phase1Idx;
137 static_assert(numPhases > 1,
138 "chemical non-equlibrium only makes sense for multiple phases");
139 static_assert(numPhases == numComponents,
140 "currently chemical non-equilibrium is only available when numPhases equals numComponents");
141 static_assert(ModelTraits::useMoles(),
142 "chemical nonequilibrium can only be calculated based on mole fractions not mass fractions");
145 using ParentType::ParentType;
154 const SubControlVolume& scv,
155 const VolumeVariables& volVars)
const
157 NumEqVector storage(0.0);
160 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
162 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
164 auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
165 storage[eqIdx] += volVars.porosity()
166 * volVars.saturation(phaseIdx)
167 * volVars.molarDensity(phaseIdx)
168 * volVars.moleFraction(phaseIdx, compIdx);
171 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
174 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
189 const Element& element,
190 const FVElementGeometry& fvGeometry,
191 const ElementVolumeVariables& elemVolVars,
192 const SubControlVolumeFace& scvf,
193 const ElementFluxVariablesCache& elemFluxVarsCache)
const
195 FluxVariables fluxVars;
196 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
197 NumEqVector flux(0.0);
200 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
202 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
203 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
206 const auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
209 const auto upwindTerm = [phaseIdx, compIdx] (
const auto& volVars)
210 {
return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
212 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
215 if (compIdx == phaseIdx)
219 flux[eqIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
221 flux[eqIdx] += diffusiveFluxes[compIdx];
224 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
228 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
243 const Element& element,
244 const FVElementGeometry& fvGeometry,
245 const ElementVolumeVariables& elemVolVars,
246 const SubControlVolume &scv)
const
248 NumEqVector source(0.0);
252 const auto& volVars = elemVolVars[scv];
253 std::array<std::array<Scalar, numComponents>, numPhases> componentIntoPhaseMassTransfer = {{{0.0},{0.0}}};
256 const Scalar characteristicLength = volVars.characteristicLength() ;
257 const Scalar factorMassTransfer = volVars.factorMassTransfer() ;
259 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
261 for(
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
263 const Scalar sherwoodNumber = volVars.sherwoodNumber(phaseIdx);
265 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
267 if (compIdx <= numPhases)
270 if (phaseIdx == compIdx)
273 const Scalar xNonEquil = volVars.moleFraction(phaseIdx, compIdx);
276 const Scalar xEquil = volVars.xEquil(phaseIdx, compIdx);
278 const Scalar diffCoeff = volVars.diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
281 const Scalar compFluxIntoOtherPhase = factorMassTransfer * (xEquil-xNonEquil)/characteristicLength * awn * volVars.molarDensity(phaseIdx) * diffCoeff * sherwoodNumber;
283 componentIntoPhaseMassTransfer[phaseIdx][compIdx] += compFluxIntoOtherPhase;
284 componentIntoPhaseMassTransfer[compIdx][compIdx] -= compFluxIntoOtherPhase;
291 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
293 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
295 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
296 source[eqIdx] += componentIntoPhaseMassTransfer[phaseIdx][compIdx];
299 if (!isfinite(source[eqIdx]))
304 if constexpr (ModelTraits::enableThermalNonEquilibrium())
308 EnergyLocalResidual::computeSourceEnergy(source,
316 source += problem.source(element, fvGeometry, elemVolVars, scv);
319 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
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
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Definition: porousmediumflow/nonequilibrium/localresidual.hh:37
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:73
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:242
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:153
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:188
Declares all properties used in Dumux.
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...