3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/3pwateroil/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_3P2CNI_LOCAL_RESIDUAL_HH
27#define DUMUX_3P2CNI_LOCAL_RESIDUAL_HH
28
32
33namespace Dumux {
39template<class TypeTag>
40class ThreePWaterOilLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
41{
42protected:
47 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
48 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using Element = typename GridView::template Codim<0>::Entity;
59
60 enum {
63
64 conti0EqIdx = Indices::conti0EqIdx,
66
67 wPhaseIdx = FluidSystem::wPhaseIdx,
68 nPhaseIdx = FluidSystem::nPhaseIdx,
69 gPhaseIdx = FluidSystem::gPhaseIdx,
70
71 wCompIdx = FluidSystem::wCompIdx,
72 nCompIdx = FluidSystem::nCompIdx,
73 };
74
76 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
77public:
78 using ParentType::ParentType;
79
92 const SubControlVolume& scv,
93 const VolumeVariables& volVars) const
94 {
95 NumEqVector storage(0.0);
96
97 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
98 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
99
100 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
101 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
102
103 // compute storage term of all components within all phases
104 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
105 {
106 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
107 {
108 int eqIdx = (compIdx == wCompIdx) ? conti0EqIdx : conti1EqIdx;
109 storage[eqIdx] += volVars.porosity()
110 * volVars.saturation(phaseIdx)
111 * massOrMoleDensity(volVars, phaseIdx)
112 * massOrMoleFraction(volVars, phaseIdx, compIdx);
113 }
114
116 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
117 }
118
120 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
121
122 return storage;
123 }
124
137 const Element& element,
138 const FVElementGeometry& fvGeometry,
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf,
141 const ElementFluxVariablesCache& elemFluxVarsCache) const
142 {
143 FluxVariables fluxVars;
144 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
145
146 // get upwind weights into local scope
147 NumEqVector flux(0.0);
148
149 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
150 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
151
152 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
153 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
154
155 // advective fluxes
156 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
157 {
158 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
159 {
160 const auto upwindTerm = [&massOrMoleDensity, &massOrMoleFraction, phaseIdx, compIdx] (const auto& volVars)
161 { return massOrMoleDensity(volVars, phaseIdx)*massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
162
163 // get equation index
164 auto eqIdx = conti0EqIdx + compIdx;
165 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
166 }
167
169 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
170 }
171
173 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
174
175 // diffusive fluxes
176 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
177 const auto diffusionFluxesWPhase = fluxVars.molecularDiffusionFlux(wPhaseIdx);
178 Scalar jNW = diffusionFluxesWPhase[nCompIdx];
179 Scalar jWW = -jNW;
180 // check for the reference system and adapt units of the diffusive flux accordingly.
181 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
182 {
183 jNW /= FluidSystem::molarMass(nCompIdx);
184 jWW /= FluidSystem::molarMass(wCompIdx);
185 }
186
187 const auto diffusionFluxesGPhase = fluxVars.molecularDiffusionFlux(gPhaseIdx);
188 Scalar jWG = diffusionFluxesGPhase[wCompIdx];
189 Scalar jNG = -jWG;
190 // check for the reference system and adapt units of the diffusive flux accordingly.
191 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
192 {
193 jWG /= FluidSystem::molarMass(wCompIdx);
194 jNG /= FluidSystem::molarMass(nCompIdx);
195 }
196
197 const auto diffusionFluxesNPhase = fluxVars.molecularDiffusionFlux(nPhaseIdx);
198 Scalar jWN = diffusionFluxesNPhase[wCompIdx];
199 Scalar jNN = -jWN;
200 // check for the reference system and adapt units of the diffusive flux accordingly.
201 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
202 {
203 jWN /= FluidSystem::molarMass(wCompIdx);
204 jNN /= FluidSystem::molarMass(nCompIdx);
205 }
206
207 flux[conti0EqIdx] += jWW+jWG+jWN;
208 flux[conti1EqIdx] += jNW+jNG+jNN;
209
210 return flux;
211 }
212};
213
214} // end namespace Dumux
215
216#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 the ThreePWaterOil fully implicit m...
Definition: porousmediumflow/3pwateroil/localresidual.hh:41
typename GetPropType< TypeTag, Properties::GridGeometry >::LocalView FVElementGeometry
Definition: porousmediumflow/3pwateroil/localresidual.hh:46
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/3pwateroil/localresidual.hh:136
@ numComponents
Definition: porousmediumflow/3pwateroil/localresidual.hh:62
@ conti0EqIdx
Index of the mass conservation equation for the water component.
Definition: porousmediumflow/3pwateroil/localresidual.hh:64
@ nPhaseIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:68
@ gPhaseIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:69
@ conti1EqIdx
Index of the mass conservation equation for the contaminant component.
Definition: porousmediumflow/3pwateroil/localresidual.hh:65
@ wPhaseIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:67
@ nCompIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:72
@ numPhases
Definition: porousmediumflow/3pwateroil/localresidual.hh:61
@ wCompIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:71
typename GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView ElementFluxVariablesCache
Definition: porousmediumflow/3pwateroil/localresidual.hh:51
GetPropType< TypeTag, Properties::VolumeVariables > VolumeVariables
Definition: porousmediumflow/3pwateroil/localresidual.hh:56
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/3pwateroil/localresidual.hh:91
typename GetPropType< TypeTag, Properties::GridGeometry >::GridView GridView
Definition: porousmediumflow/3pwateroil/localresidual.hh:53
typename GridView::template Codim< 0 >::Entity Element
Definition: porousmediumflow/3pwateroil/localresidual.hh:54
GetPropType< TypeTag, Properties::Problem > Problem
Definition: porousmediumflow/3pwateroil/localresidual.hh:44
typename GetPropType< TypeTag, Properties::ModelTraits >::Indices Indices
Definition: porousmediumflow/3pwateroil/localresidual.hh:52
typename GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView ElementVolumeVariables
Definition: porousmediumflow/3pwateroil/localresidual.hh:55
GetPropType< TypeTag, Properties::FluxVariables > FluxVariables
Definition: porousmediumflow/3pwateroil/localresidual.hh:50
GetPropType< TypeTag, Properties::BaseLocalResidual > ParentType
Definition: porousmediumflow/3pwateroil/localresidual.hh:43
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: porousmediumflow/3pwateroil/localresidual.hh:58
typename FVElementGeometry::SubControlVolumeFace SubControlVolumeFace
Definition: porousmediumflow/3pwateroil/localresidual.hh:48
GetPropType< TypeTag, Properties::EnergyLocalResidual > EnergyLocalResidual
Definition: porousmediumflow/3pwateroil/localresidual.hh:57
Dumux::NumEqVector< GetPropType< TypeTag, Properties::PrimaryVariables > > NumEqVector
Definition: porousmediumflow/3pwateroil/localresidual.hh:49
static constexpr bool useMoles
Property that defines whether mole or mass fractions are used.
Definition: porousmediumflow/3pwateroil/localresidual.hh:76
typename FVElementGeometry::SubControlVolume SubControlVolume
Definition: porousmediumflow/3pwateroil/localresidual.hh:47
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: porousmediumflow/3pwateroil/localresidual.hh:45
Declares all properties used in Dumux.