26#ifndef DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
27#define DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
30#include <dune/common/exceptions.hh>
41template<
class TypeTag>
49 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
50 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
55 using Element =
typename GridView::template Codim<0>::Entity;
61 using Indices =
typename ModelTraits::Indices;
63 static constexpr int numPhases = ModelTraits::numFluidPhases();
64 static constexpr int numComponents = ModelTraits::numFluidComponents();
65 static constexpr bool useMoles = ModelTraits::useMoles();
67 enum { conti0EqIdx = Indices::conti0EqIdx };
70 static constexpr int replaceCompEqIdx = ModelTraits::replaceCompEqIdx();
71 static constexpr bool useTotalMoleOrMassBalance = replaceCompEqIdx < numComponents;
74 using ParentType::ParentType;
88 const SubControlVolume& scv,
89 const VolumeVariables& volVars)
const
91 NumEqVector storage(0.0);
93 const auto massOrMoleDensity = [](
const auto& volVars,
const int phaseIdx)
94 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
96 const auto massOrMoleFraction= [](
const auto& volVars,
const int phaseIdx,
const int compIdx)
97 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
100 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
102 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
104 auto eqIdx = conti0EqIdx + compIdx;
105 if (eqIdx != replaceCompEqIdx)
106 storage[eqIdx] += volVars.porosity()
107 * volVars.saturation(phaseIdx)
108 * massOrMoleDensity(volVars, phaseIdx)
113 if (useTotalMoleOrMassBalance)
114 storage[replaceCompEqIdx] += massOrMoleDensity(volVars, phaseIdx)
116 * volVars.saturation(phaseIdx);
119 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
123 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
140 const Element& element,
141 const FVElementGeometry& fvGeometry,
142 const ElementVolumeVariables& elemVolVars,
143 const SubControlVolumeFace& scvf,
144 const ElementFluxVariablesCache& elemFluxVarsCache)
const
146 FluxVariables fluxVars;
147 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
148 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
150 NumEqVector flux(0.0);
152 const auto massOrMoleDensity = [](
const auto& volVars,
const int phaseIdx)
153 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
155 const auto massOrMoleFraction= [](
const auto& volVars,
const int phaseIdx,
const int compIdx)
156 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
159 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
161 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
162 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
165 const auto eqIdx = conti0EqIdx + compIdx;
168 const auto upwindTerm = [&massOrMoleDensity, &
massOrMoleFraction, phaseIdx, compIdx] (
const auto& volVars)
169 {
return massOrMoleDensity(volVars, phaseIdx)*
massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
171 if (eqIdx != replaceCompEqIdx)
172 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
175 if(eqIdx != replaceCompEqIdx)
179 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx)
180 : diffusiveFluxes[compIdx];
182 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]
183 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
185 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
190 if (useTotalMoleOrMassBalance)
193 const auto upwindTerm = [&massOrMoleDensity, phaseIdx] (
const auto& volVars)
194 {
return massOrMoleDensity(volVars, phaseIdx)*volVars.mobility(phaseIdx); };
196 flux[replaceCompEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
198 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
202 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
204 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]
205 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
207 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
212 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
216 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
223 {
return static_cast<Implementation *
> (
this); }
226 {
return static_cast<const Implementation *
> (
this); }
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 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
Element-wise calculation of the local residual for problems using compositional fully implicit model.
Definition: porousmediumflow/compositional/localresidual.hh:43
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:139
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:87
Implementation * asImp_()
Definition: porousmediumflow/compositional/localresidual.hh:222
const Implementation * asImp_() const
Definition: porousmediumflow/compositional/localresidual.hh:225
Declares all properties used in Dumux.