version 3.10-dev
porousmediumflow/compositional/localresidual.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
15#define DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
16
17#include <vector>
18#include <dune/common/exceptions.hh>
23
24namespace Dumux {
25
31template<class TypeTag>
32class CompositionalLocalResidual: public GetPropType<TypeTag, Properties::BaseLocalResidual>
33{
38 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
39 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
40 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
43 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
45 using Element = typename GridView::template Codim<0>::Entity;
46 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
51 using Indices = typename ModelTraits::Indices;
52
53 static constexpr int numPhases = ModelTraits::numFluidPhases();
54 static constexpr int numComponents = ModelTraits::numFluidComponents();
55 static constexpr bool useMoles = ModelTraits::useMoles();
56
57 enum { conti0EqIdx = Indices::conti0EqIdx };
58
60 static constexpr int replaceCompEqIdx = ModelTraits::replaceCompEqIdx();
61 static constexpr bool useTotalMoleOrMassBalance = replaceCompEqIdx < numComponents;
62
63public:
64 using ParentType::ParentType;
65
77 NumEqVector computeStorage(const Problem& problem,
78 const SubControlVolume& scv,
79 const VolumeVariables& volVars) const
80 {
81 NumEqVector storage(0.0);
82
83 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
84 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
85
86 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
87 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
88
89 // compute storage term of all components within all phases
90 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
91 {
92 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
93 {
94 auto eqIdx = conti0EqIdx + compIdx;
95 if (eqIdx != replaceCompEqIdx)
96 storage[eqIdx] += volVars.porosity()
97 * volVars.saturation(phaseIdx)
98 * massOrMoleDensity(volVars, phaseIdx)
99 * massOrMoleFraction(volVars, phaseIdx, compIdx);
100 }
101
102 // in case one balance is substituted by the total mole balance
103 if (useTotalMoleOrMassBalance)
104 storage[replaceCompEqIdx] += massOrMoleDensity(volVars, phaseIdx)
105 * volVars.porosity()
106 * volVars.saturation(phaseIdx);
107
109 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx);
110 }
111
113 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
114
115 return storage;
116 }
117
129 NumEqVector computeFlux(const Problem& problem,
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
134 const ElementFluxVariablesCache& elemFluxVarsCache) const
135 {
136 FluxVariables fluxVars;
137 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
138 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
139 // get upwind weights into local scope
140 NumEqVector flux(0.0);
141
142 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
143 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
144
145 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
146 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
147
148 // advective fluxes
149 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
150 {
151 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
152 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
153 {
154 // get equation index
155 const auto eqIdx = conti0EqIdx + compIdx;
156
157 // the physical quantities for which we perform upwinding
158 const auto upwindTerm = [&massOrMoleDensity, &massOrMoleFraction, phaseIdx, compIdx] (const auto& volVars)
159 { return massOrMoleDensity(volVars, phaseIdx)*massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
160
161 if (eqIdx != replaceCompEqIdx)
162 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
163
164 // diffusive fluxes (only for the component balances)
165 if(eqIdx != replaceCompEqIdx)
166 {
167 //check for the reference system and adapt units of the diffusive flux accordingly.
168 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
169 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx)
170 : diffusiveFluxes[compIdx];
171 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
172 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]
173 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
174 else
175 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
176 }
177 }
178
179 // in case one balance is substituted by the total mole balance
180 if (useTotalMoleOrMassBalance)
181 {
182 // the physical quantities for which we perform upwinding
183 const auto upwindTerm = [&massOrMoleDensity, phaseIdx] (const auto& volVars)
184 { return massOrMoleDensity(volVars, phaseIdx)*volVars.mobility(phaseIdx); };
185
186 flux[replaceCompEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
187
188 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
189 {
190 //check for the reference system and adapt units of the diffusive flux accordingly.
191 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
192 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
193 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
194 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]
195 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
196 else
197 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
198 }
199 }
200
202 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
203
204 if constexpr (ModelTraits::enableCompositionalDispersion())
205 {
206 if constexpr (FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::box && numPhases == 1)
207 {
208 const auto dispersionFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
209 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
210 {
211 flux[compIdx] += dispersionFluxes[compIdx];
212 }
213 }
214 else
215 DUNE_THROW(Dune::NotImplemented, "Dispersion Fluxes are only implemented for single phase flows using the Box method.");
216 }
217
218 }
219
221 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
222 EnergyLocalResidual::heatDispersionFlux(flux, fluxVars);
223
224 return flux;
225 }
226
227protected:
228 Implementation *asImp_()
229 { return static_cast<Implementation *> (this); }
230
231 const Implementation *asImp_() const
232 { return static_cast<const Implementation *> (this); }
233};
234
235} // end namespace Dumux
236
237#endif
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
Definition: adapt.hh:17
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...