14#ifndef DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
15#define DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
18#include <dune/common/exceptions.hh>
31template<
class TypeTag>
39 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
40 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
45 using Element =
typename GridView::template Codim<0>::Entity;
51 using Indices =
typename ModelTraits::Indices;
53 static constexpr int numPhases = ModelTraits::numFluidPhases();
54 static constexpr int numComponents = ModelTraits::numFluidComponents();
55 static constexpr bool useMoles = ModelTraits::useMoles();
57 enum { conti0EqIdx = Indices::conti0EqIdx };
60 static constexpr int replaceCompEqIdx = ModelTraits::replaceCompEqIdx();
61 static constexpr bool useTotalMoleOrMassBalance = replaceCompEqIdx < numComponents;
64 using ParentType::ParentType;
78 const SubControlVolume& scv,
79 const VolumeVariables& volVars)
const
81 NumEqVector storage(0.0);
83 const auto massOrMoleDensity = [](
const auto& volVars,
const int phaseIdx)
84 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
86 const auto massOrMoleFraction= [](
const auto& volVars,
const int phaseIdx,
const int compIdx)
87 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
90 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
92 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
94 auto eqIdx = conti0EqIdx + compIdx;
95 if (eqIdx != replaceCompEqIdx)
96 storage[eqIdx] += volVars.porosity()
97 * volVars.saturation(phaseIdx)
98 * massOrMoleDensity(volVars, phaseIdx)
103 if (useTotalMoleOrMassBalance)
104 storage[replaceCompEqIdx] += massOrMoleDensity(volVars, phaseIdx)
106 * volVars.saturation(phaseIdx);
109 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
113 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
134 const ElementFluxVariablesCache& elemFluxVarsCache)
const
136 FluxVariables fluxVars;
137 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
138 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
140 NumEqVector flux(0.0);
142 const auto massOrMoleDensity = [](
const auto& volVars,
const int phaseIdx)
143 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
145 const auto massOrMoleFraction= [](
const auto& volVars,
const int phaseIdx,
const int compIdx)
146 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
149 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
151 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
152 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
155 const auto eqIdx = conti0EqIdx + compIdx;
158 const auto upwindTerm = [&massOrMoleDensity, &
massOrMoleFraction, phaseIdx, compIdx] (
const auto& volVars)
159 {
return massOrMoleDensity(volVars, phaseIdx)*
massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
161 if (eqIdx != replaceCompEqIdx)
162 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
165 if(eqIdx != replaceCompEqIdx)
169 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx)
170 : diffusiveFluxes[compIdx];
172 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]
173 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
175 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
180 if (useTotalMoleOrMassBalance)
183 const auto upwindTerm = [&massOrMoleDensity, phaseIdx] (
const auto& volVars)
184 {
return massOrMoleDensity(volVars, phaseIdx)*volVars.mobility(phaseIdx); };
186 flux[replaceCompEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
188 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
192 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
194 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]
195 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
197 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
202 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
204 if constexpr (ModelTraits::enableCompositionalDispersion())
208 const auto dispersionFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
209 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
211 flux[compIdx] += dispersionFluxes[compIdx];
215 DUNE_THROW(Dune::NotImplemented,
"Dispersion Fluxes are only implemented for single phase flows using the Box method.");
221 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
222 EnergyLocalResidual::heatDispersionFlux(flux, fluxVars);
229 {
return static_cast<Implementation *
> (
this); }
232 {
return static_cast<const Implementation *
> (
this); }
Element-wise calculation of the local residual for problems using compositional fully implicit model.
Definition: porousmediumflow/compositional/localresidual.hh:33
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluates the total flux of all conservation quantities over a face of a sub-control volume.
Definition: porousmediumflow/compositional/localresidual.hh:129
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluates the amount of all conservation quantities (e.g. phase mass) within a sub-control volume.
Definition: porousmediumflow/compositional/localresidual.hh:77
Implementation * asImp_()
Definition: porousmediumflow/compositional/localresidual.hh:228
const Implementation * asImp_() const
Definition: porousmediumflow/compositional/localresidual.hh:231
Defines all properties used in Dumux.
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:34
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:54
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147
A helper to deduce a vector with the same size as numbers of equations.