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