version 3.10-dev
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// 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_3P2CNI_LOCAL_RESIDUAL_HH
15#define DUMUX_3P2CNI_LOCAL_RESIDUAL_HH
16
20
21namespace Dumux {
27template<class TypeTag>
28class ThreePWaterOilLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
29{
30protected:
35 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
36 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
42 using Element = typename GridView::template Codim<0>::Entity;
47
48 enum {
51
52 conti0EqIdx = Indices::conti0EqIdx,
54
55 wPhaseIdx = FluidSystem::wPhaseIdx,
56 nPhaseIdx = FluidSystem::nPhaseIdx,
57 gPhaseIdx = FluidSystem::gPhaseIdx,
58
59 wCompIdx = FluidSystem::wCompIdx,
60 nCompIdx = FluidSystem::nCompIdx,
61 };
62
64 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
65public:
66 using ParentType::ParentType;
67
80 const SubControlVolume& scv,
81 const VolumeVariables& volVars) const
82 {
83 NumEqVector storage(0.0);
84
85 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
86 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
87
88 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
89 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
90
91 // compute storage term of all components within all phases
92 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
93 {
94 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
95 {
96 int eqIdx = (compIdx == wCompIdx) ? conti0EqIdx : conti1EqIdx;
97 storage[eqIdx] += volVars.porosity()
98 * volVars.saturation(phaseIdx)
99 * massOrMoleDensity(volVars, phaseIdx)
100 * massOrMoleFraction(volVars, phaseIdx, compIdx);
101 }
102
104 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx);
105 }
106
108 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
109
110 return storage;
111 }
112
125 const Element& element,
126 const FVElementGeometry& fvGeometry,
127 const ElementVolumeVariables& elemVolVars,
128 const SubControlVolumeFace& scvf,
129 const ElementFluxVariablesCache& elemFluxVarsCache) const
130 {
131 FluxVariables fluxVars;
132 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
133
134 // get upwind weights into local scope
135 NumEqVector flux(0.0);
136
137 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
138 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
139
140 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
141 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
142
143 // advective fluxes
144 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
145 {
146 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
147 {
148 const auto upwindTerm = [&massOrMoleDensity, &massOrMoleFraction, phaseIdx, compIdx] (const auto& volVars)
149 { return massOrMoleDensity(volVars, phaseIdx)*massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
150
151 // get equation index
152 auto eqIdx = conti0EqIdx + compIdx;
153 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
154 }
155
157 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
158 }
159
161 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
162
163 // diffusive fluxes
164 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
165 const auto diffusionFluxesWPhase = fluxVars.molecularDiffusionFlux(wPhaseIdx);
166 Scalar jNW = diffusionFluxesWPhase[nCompIdx];
167 Scalar jWW = -jNW;
168 // check for the reference system and adapt units of the diffusive flux accordingly.
169 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
170 {
171 jNW /= FluidSystem::molarMass(nCompIdx);
172 jWW /= FluidSystem::molarMass(wCompIdx);
173 }
174
175 const auto diffusionFluxesGPhase = fluxVars.molecularDiffusionFlux(gPhaseIdx);
176 Scalar jWG = diffusionFluxesGPhase[wCompIdx];
177 Scalar jNG = -jWG;
178 // check for the reference system and adapt units of the diffusive flux accordingly.
179 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
180 {
181 jWG /= FluidSystem::molarMass(wCompIdx);
182 jNG /= FluidSystem::molarMass(nCompIdx);
183 }
184
185 const auto diffusionFluxesNPhase = fluxVars.molecularDiffusionFlux(nPhaseIdx);
186 Scalar jWN = diffusionFluxesNPhase[wCompIdx];
187 Scalar jNN = -jWN;
188 // check for the reference system and adapt units of the diffusive flux accordingly.
189 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
190 {
191 jWN /= FluidSystem::molarMass(wCompIdx);
192 jNN /= FluidSystem::molarMass(nCompIdx);
193 }
194
195 flux[conti0EqIdx] += jWW+jWG+jWN;
196 flux[conti1EqIdx] += jNW+jNG+jNN;
197
198 return flux;
199 }
200};
201
202} // end namespace Dumux
203
204#endif
Element-wise calculation of the local residual for problems using the ThreePWaterOil fully implicit m...
Definition: porousmediumflow/3pwateroil/localresidual.hh:29
typename GetPropType< TypeTag, Properties::GridGeometry >::LocalView FVElementGeometry
Definition: porousmediumflow/3pwateroil/localresidual.hh:34
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:124
typename GetPropType< TypeTag, Properties::GridFluxVariablesCache >::LocalView ElementFluxVariablesCache
Definition: porousmediumflow/3pwateroil/localresidual.hh:39
GetPropType< TypeTag, Properties::VolumeVariables > VolumeVariables
Definition: porousmediumflow/3pwateroil/localresidual.hh:44
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:79
typename GetPropType< TypeTag, Properties::GridGeometry >::GridView GridView
Definition: porousmediumflow/3pwateroil/localresidual.hh:41
typename GridView::template Codim< 0 >::Entity Element
Definition: porousmediumflow/3pwateroil/localresidual.hh:42
GetPropType< TypeTag, Properties::Problem > Problem
Definition: porousmediumflow/3pwateroil/localresidual.hh:32
typename GetPropType< TypeTag, Properties::ModelTraits >::Indices Indices
Definition: porousmediumflow/3pwateroil/localresidual.hh:40
typename GetPropType< TypeTag, Properties::GridVolumeVariables >::LocalView ElementVolumeVariables
Definition: porousmediumflow/3pwateroil/localresidual.hh:43
GetPropType< TypeTag, Properties::FluxVariables > FluxVariables
Definition: porousmediumflow/3pwateroil/localresidual.hh:38
@ numComponents
Definition: porousmediumflow/3pwateroil/localresidual.hh:50
@ conti0EqIdx
Index of the mass conservation equation for the water component.
Definition: porousmediumflow/3pwateroil/localresidual.hh:52
@ nPhaseIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:56
@ gPhaseIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:57
@ conti1EqIdx
Index of the mass conservation equation for the contaminant component.
Definition: porousmediumflow/3pwateroil/localresidual.hh:53
@ wPhaseIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:55
@ nCompIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:60
@ numPhases
Definition: porousmediumflow/3pwateroil/localresidual.hh:49
@ wCompIdx
Definition: porousmediumflow/3pwateroil/localresidual.hh:59
GetPropType< TypeTag, Properties::BaseLocalResidual > ParentType
Definition: porousmediumflow/3pwateroil/localresidual.hh:31
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: porousmediumflow/3pwateroil/localresidual.hh:46
typename FVElementGeometry::SubControlVolumeFace SubControlVolumeFace
Definition: porousmediumflow/3pwateroil/localresidual.hh:36
GetPropType< TypeTag, Properties::EnergyLocalResidual > EnergyLocalResidual
Definition: porousmediumflow/3pwateroil/localresidual.hh:45
Dumux::NumEqVector< GetPropType< TypeTag, Properties::PrimaryVariables > > NumEqVector
Definition: porousmediumflow/3pwateroil/localresidual.hh:37
static constexpr bool useMoles
Property that defines whether mole or mass fractions are used.
Definition: porousmediumflow/3pwateroil/localresidual.hh:64
typename FVElementGeometry::SubControlVolume SubControlVolume
Definition: porousmediumflow/3pwateroil/localresidual.hh:35
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: porousmediumflow/3pwateroil/localresidual.hh:33
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
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...