26#ifndef DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
27#define DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
30#include <dune/common/exceptions.hh>
42template<
class TypeTag>
50 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
51 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
56 using Element =
typename GridView::template Codim<0>::Entity;
62 using Indices =
typename ModelTraits::Indices;
64 static constexpr int numPhases = ModelTraits::numFluidPhases();
65 static constexpr int numComponents = ModelTraits::numFluidComponents();
66 static constexpr bool useMoles = ModelTraits::useMoles();
68 enum { conti0EqIdx = Indices::conti0EqIdx };
71 static constexpr int replaceCompEqIdx = ModelTraits::replaceCompEqIdx();
72 static constexpr bool useTotalMoleOrMassBalance = replaceCompEqIdx < numComponents;
75 using ParentType::ParentType;
89 const SubControlVolume& scv,
90 const VolumeVariables& volVars)
const
92 NumEqVector storage(0.0);
94 const auto massOrMoleDensity = [](
const auto& volVars,
const int phaseIdx)
95 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
97 const auto massOrMoleFraction= [](
const auto& volVars,
const int phaseIdx,
const int compIdx)
98 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
101 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
103 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
105 auto eqIdx = conti0EqIdx + compIdx;
106 if (eqIdx != replaceCompEqIdx)
107 storage[eqIdx] += volVars.porosity()
108 * volVars.saturation(phaseIdx)
109 * massOrMoleDensity(volVars, phaseIdx)
114 if (useTotalMoleOrMassBalance)
115 storage[replaceCompEqIdx] += massOrMoleDensity(volVars, phaseIdx)
117 * volVars.saturation(phaseIdx);
120 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
124 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
141 const Element& element,
142 const FVElementGeometry& fvGeometry,
143 const ElementVolumeVariables& elemVolVars,
144 const SubControlVolumeFace& scvf,
145 const ElementFluxVariablesCache& elemFluxVarsCache)
const
147 FluxVariables fluxVars;
148 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
149 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
151 NumEqVector flux(0.0);
153 const auto massOrMoleDensity = [](
const auto& volVars,
const int phaseIdx)
154 {
return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
156 const auto massOrMoleFraction= [](
const auto& volVars,
const int phaseIdx,
const int compIdx)
157 {
return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
160 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
162 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
163 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
166 const auto eqIdx = conti0EqIdx + compIdx;
169 const auto upwindTerm = [&massOrMoleDensity, &
massOrMoleFraction, phaseIdx, compIdx] (
const auto& volVars)
170 {
return massOrMoleDensity(volVars, phaseIdx)*
massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
172 if (eqIdx != replaceCompEqIdx)
173 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
176 if(eqIdx != replaceCompEqIdx)
180 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx)
181 : diffusiveFluxes[compIdx];
183 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]
184 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
186 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
191 if (useTotalMoleOrMassBalance)
194 const auto upwindTerm = [&massOrMoleDensity, phaseIdx] (
const auto& volVars)
195 {
return massOrMoleDensity(volVars, phaseIdx)*volVars.mobility(phaseIdx); };
197 flux[replaceCompEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
199 for(
int compIdx = 0; compIdx < numComponents; ++compIdx)
203 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
205 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]
206 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
208 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
213 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
217 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
224 {
return static_cast<Implementation *
> (
this); }
227 {
return static_cast<const Implementation *
> (
this); }
A helper to deduce a vector with the same size as numbers of equations.
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
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Element-wise calculation of the local residual for problems using compositional fully implicit model.
Definition: porousmediumflow/compositional/localresidual.hh:44
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:140
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:88
Implementation * asImp_()
Definition: porousmediumflow/compositional/localresidual.hh:223
const Implementation * asImp_() const
Definition: porousmediumflow/compositional/localresidual.hh:226
Declares all properties used in Dumux.