26#ifndef DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
27#define DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
30#include <dune/common/exceptions.hh>
43template<
class TypeTag>
51 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
52 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
57 using Element =
typename GridView::template Codim<0>::Entity;
63 using Indices =
typename ModelTraits::Indices;
65 static constexpr int numPhases = ModelTraits::numFluidPhases();
66 static constexpr int numComponents = ModelTraits::numFluidComponents();
67 static constexpr bool useMoles = ModelTraits::useMoles();
69 enum { conti0EqIdx = Indices::conti0EqIdx };
72 static constexpr int replaceCompEqIdx = ModelTraits::replaceCompEqIdx();
73 static constexpr bool useTotalMoleOrMassBalance = replaceCompEqIdx < numComponents;
76 using ParentType::ParentType;
90 const SubControlVolume& scv,
91 const VolumeVariables& volVars)
const
93 NumEqVector storage(0.0);
95 const auto massOrMoleDensity = [](
const auto& volVars,
const int phaseIdx)
96 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
98 const auto massOrMoleFraction= [](
const auto& volVars,
const int phaseIdx,
const int compIdx)
99 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
102 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
104 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
106 auto eqIdx = conti0EqIdx + compIdx;
107 if (eqIdx != replaceCompEqIdx)
108 storage[eqIdx] += volVars.porosity()
109 * volVars.saturation(phaseIdx)
110 * massOrMoleDensity(volVars, phaseIdx)
115 if (useTotalMoleOrMassBalance)
116 storage[replaceCompEqIdx] += massOrMoleDensity(volVars, phaseIdx)
118 * volVars.saturation(phaseIdx);
121 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
125 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
142 const Element& element,
143 const FVElementGeometry& fvGeometry,
144 const ElementVolumeVariables& elemVolVars,
145 const SubControlVolumeFace& scvf,
146 const ElementFluxVariablesCache& elemFluxVarsCache)
const
148 FluxVariables fluxVars;
149 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
150 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
152 NumEqVector flux(0.0);
154 const auto massOrMoleDensity = [](
const auto& volVars,
const int phaseIdx)
155 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
157 const auto massOrMoleFraction= [](
const auto& volVars,
const int phaseIdx,
const int compIdx)
158 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
161 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
163 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
164 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
167 const auto eqIdx = conti0EqIdx + compIdx;
170 const auto upwindTerm = [&massOrMoleDensity, &
massOrMoleFraction, phaseIdx, compIdx] (
const auto& volVars)
171 {
return massOrMoleDensity(volVars, phaseIdx)*
massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
173 if (eqIdx != replaceCompEqIdx)
174 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
177 if(eqIdx != replaceCompEqIdx)
181 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx)
182 : diffusiveFluxes[compIdx];
184 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]
185 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
187 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
192 if (useTotalMoleOrMassBalance)
195 const auto upwindTerm = [&massOrMoleDensity, phaseIdx] (
const auto& volVars)
196 {
return massOrMoleDensity(volVars, phaseIdx)*volVars.mobility(phaseIdx); };
198 flux[replaceCompEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
200 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
204 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
206 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]
207 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
209 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
214 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
216 if constexpr (ModelTraits::enableCompositionalDispersion())
220 const auto dispersionFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
221 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
223 flux[compIdx] += dispersionFluxes[compIdx];
227 DUNE_THROW(Dune::NotImplemented,
"Dispersion Fluxes are only implemented for single phase flows using the Box method.");
233 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
234 EnergyLocalResidual::heatDispersionFlux(flux, fluxVars);
241 {
return static_cast<Implementation *
> (
this); }
244 {
return static_cast<const Implementation *
> (
this); }
A helper to deduce a vector with the same size as numbers of equations.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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:66
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
Element-wise calculation of the local residual for problems using compositional fully implicit model.
Definition: porousmediumflow/compositional/localresidual.hh:45
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:141
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:89
Implementation * asImp_()
Definition: porousmediumflow/compositional/localresidual.hh:240
const Implementation * asImp_() const
Definition: porousmediumflow/compositional/localresidual.hh:243
Declares all properties used in Dumux.