3.6-git
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>
35
36namespace Dumux {
37
43template<class TypeTag>
44class CompositionalLocalResidual: public GetPropType<TypeTag, Properties::BaseLocalResidual>
45{
50 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
51 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
52 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
55 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
57 using Element = typename GridView::template Codim<0>::Entity;
58 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
63 using Indices = typename ModelTraits::Indices;
64
65 static constexpr int numPhases = ModelTraits::numFluidPhases();
66 static constexpr int numComponents = ModelTraits::numFluidComponents();
67 static constexpr bool useMoles = ModelTraits::useMoles();
68
69 enum { conti0EqIdx = Indices::conti0EqIdx };
70
72 static constexpr int replaceCompEqIdx = ModelTraits::replaceCompEqIdx();
73 static constexpr bool useTotalMoleOrMassBalance = replaceCompEqIdx < numComponents;
74
75public:
76 using ParentType::ParentType;
77
89 NumEqVector computeStorage(const Problem& problem,
90 const SubControlVolume& scv,
91 const VolumeVariables& volVars) const
92 {
93 NumEqVector storage(0.0);
94
95 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
96 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
97
98 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
99 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
100
101 // compute storage term of all components within all phases
102 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
103 {
104 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
105 {
106 auto eqIdx = conti0EqIdx + compIdx;
107 if (eqIdx != replaceCompEqIdx)
108 storage[eqIdx] += volVars.porosity()
109 * volVars.saturation(phaseIdx)
110 * massOrMoleDensity(volVars, phaseIdx)
111 * massOrMoleFraction(volVars, phaseIdx, compIdx);
112 }
113
114 // in case one balance is substituted by the total mole balance
115 if (useTotalMoleOrMassBalance)
116 storage[replaceCompEqIdx] += massOrMoleDensity(volVars, phaseIdx)
117 * volVars.porosity()
118 * volVars.saturation(phaseIdx);
119
121 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
122 }
123
125 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
126
127 return storage;
128 }
129
141 NumEqVector computeFlux(const Problem& problem,
142 const Element& element,
143 const FVElementGeometry& fvGeometry,
144 const ElementVolumeVariables& elemVolVars,
145 const SubControlVolumeFace& scvf,
146 const ElementFluxVariablesCache& elemFluxVarsCache) const
147 {
148 FluxVariables fluxVars;
149 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
150 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
151 // get upwind weights into local scope
152 NumEqVector flux(0.0);
153
154 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
155 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
156
157 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
158 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
159
160 // advective fluxes
161 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
162 {
163 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
164 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
165 {
166 // get equation index
167 const auto eqIdx = conti0EqIdx + compIdx;
168
169 // the physical quantities for which we perform upwinding
170 const auto upwindTerm = [&massOrMoleDensity, &massOrMoleFraction, phaseIdx, compIdx] (const auto& volVars)
171 { return massOrMoleDensity(volVars, phaseIdx)*massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
172
173 if (eqIdx != replaceCompEqIdx)
174 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
175
176 // diffusive fluxes (only for the component balances)
177 if(eqIdx != replaceCompEqIdx)
178 {
179 //check for the reference system and adapt units of the diffusive flux accordingly.
180 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
181 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx)
182 : diffusiveFluxes[compIdx];
183 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
184 flux[eqIdx] += useMoles ? diffusiveFluxes[compIdx]
185 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
186 else
187 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
188 }
189 }
190
191 // in case one balance is substituted by the total mole balance
192 if (useTotalMoleOrMassBalance)
193 {
194 // the physical quantities for which we perform upwinding
195 const auto upwindTerm = [&massOrMoleDensity, phaseIdx] (const auto& volVars)
196 { return massOrMoleDensity(volVars, phaseIdx)*volVars.mobility(phaseIdx); };
197
198 flux[replaceCompEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
199
200 for(int compIdx = 0; compIdx < numComponents; ++compIdx)
201 {
202 //check for the reference system and adapt units of the diffusive flux accordingly.
203 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
204 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx) : diffusiveFluxes[compIdx];
205 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
206 flux[replaceCompEqIdx] += useMoles ? diffusiveFluxes[compIdx]
207 : diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
208 else
209 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
210 }
211 }
212
214 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
215
216 if constexpr (ModelTraits::enableCompositionalDispersion())
217 {
218 if constexpr (FVElementGeometry::GridGeometry::discMethod == DiscretizationMethods::box && numPhases == 1)
219 {
220 const auto dispersionFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
221 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
222 {
223 flux[compIdx] += dispersionFluxes[compIdx];
224 }
225 }
226 else
227 DUNE_THROW(Dune::NotImplemented, "Dispersion Fluxes are only implemented for single phase flows using the Box method.");
228 }
229
230 }
231
233 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
234 EnergyLocalResidual::heatDispersionFlux(flux, fluxVars);
235
236 return flux;
237 }
238
239protected:
240 Implementation *asImp_()
241 { return static_cast<Implementation *> (this); }
242
243 const Implementation *asImp_() const
244 { return static_cast<const Implementation *> (this); }
245};
246
247} // end namespace Dumux
248
249#endif
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.