3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
27#define DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
28
29#include <vector>
30#include <dune/common/exceptions.hh>
34
35namespace Dumux {
36
42template<class TypeTag>
43class CompositionalLocalResidual: public GetPropType<TypeTag, Properties::BaseLocalResidual>
44{
49 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
50 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
51 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
56 using Element = typename GridView::template Codim<0>::Entity;
57 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
62 using Indices = typename ModelTraits::Indices;
63
64 static constexpr int numPhases = ModelTraits::numFluidPhases();
65 static constexpr int numComponents = ModelTraits::numFluidComponents();
66 static constexpr bool useMoles = ModelTraits::useMoles();
67
68 enum { conti0EqIdx = Indices::conti0EqIdx };
69
71 static constexpr int replaceCompEqIdx = ModelTraits::replaceCompEqIdx();
72 static constexpr bool useTotalMoleOrMassBalance = replaceCompEqIdx < numComponents;
73
74public:
75 using ParentType::ParentType;
76
88 NumEqVector computeStorage(const Problem& problem,
89 const SubControlVolume& scv,
90 const VolumeVariables& volVars) const
91 {
92 NumEqVector storage(0.0);
93
94 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
95 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
96
97 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
98 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
99
100 // compute storage term of all components within all phases
101 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
102 {
103 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
104 {
105 auto eqIdx = conti0EqIdx + compIdx;
106 if (eqIdx != replaceCompEqIdx)
107 storage[eqIdx] += volVars.porosity()
108 * volVars.saturation(phaseIdx)
109 * massOrMoleDensity(volVars, phaseIdx)
110 * massOrMoleFraction(volVars, phaseIdx, compIdx);
111 }
112
113 // in case one balance is substituted by the total mole balance
114 if (useTotalMoleOrMassBalance)
115 storage[replaceCompEqIdx] += massOrMoleDensity(volVars, phaseIdx)
116 * volVars.porosity()
117 * volVars.saturation(phaseIdx);
118
120 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
121 }
122
124 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
125
126 return storage;
127 }
128
140 NumEqVector computeFlux(const Problem& problem,
141 const Element& element,
142 const FVElementGeometry& fvGeometry,
143 const ElementVolumeVariables& elemVolVars,
144 const SubControlVolumeFace& scvf,
145 const ElementFluxVariablesCache& elemFluxVarsCache) const
146 {
147 FluxVariables fluxVars;
148 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
149 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
150 // get upwind weights into local scope
151 NumEqVector flux(0.0);
152
153 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
154 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
155
156 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
157 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
158
159 // advective fluxes
160 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
161 {
162 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
163 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
164 {
165 // get equation index
166 const auto eqIdx = conti0EqIdx + compIdx;
167
168 // the physical quantities for which we perform upwinding
169 const auto upwindTerm = [&massOrMoleDensity, &massOrMoleFraction, phaseIdx, compIdx] (const auto& volVars)
170 { return massOrMoleDensity(volVars, phaseIdx)*massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
171
172 if (eqIdx != replaceCompEqIdx)
173 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
174
175 // diffusive fluxes (only for the component balances)
176 if(eqIdx != replaceCompEqIdx)
177 {
178 //check for the reference system and adapt units of the diffusive flux accordingly.
179 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
180 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx)
181 : diffusiveFluxes[compIdx];
182 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
183 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]
184 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
185 else
186 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
187 }
188 }
189
190 // in case one balance is substituted by the total mole balance
191 if (useTotalMoleOrMassBalance)
192 {
193 // the physical quantities for which we perform upwinding
194 const auto upwindTerm = [&massOrMoleDensity, phaseIdx] (const auto& volVars)
195 { return massOrMoleDensity(volVars, phaseIdx)*volVars.mobility(phaseIdx); };
196
197 flux[replaceCompEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
198
199 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
200 {
201 //check for the reference system and adapt units of the diffusive flux accordingly.
202 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
203 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
204 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
205 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]
206 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
207 else
208 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
209 }
210 }
211
213 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
214 }
215
217 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
218
219 return flux;
220 }
221
222protected:
223 Implementation *asImp_()
224 { return static_cast<Implementation *> (this); }
225
226 const Implementation *asImp_() const
227 { return static_cast<const Implementation *> (this); }
228};
229
230} // end namespace Dumux
231
232#endif
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
Definition: adapt.hh:29
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.