3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/3p3c/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_3P3C_LOCAL_RESIDUAL_HH
27#define DUMUX_3P3C_LOCAL_RESIDUAL_HH
28
32
33namespace Dumux
34{
40template<class TypeTag>
41class ThreePThreeCLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
42{
46 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
47 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
48 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
51 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
54 using Element = typename GridView::template Codim<0>::Entity;
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
59
60 enum {
63
64 contiWEqIdx = Indices::conti0EqIdx + FluidSystem::wPhaseIdx,
65 contiNEqIdx = Indices::conti0EqIdx + FluidSystem::nPhaseIdx,
66 contiGEqIdx = Indices::conti0EqIdx + FluidSystem::gPhaseIdx,
67
68 wPhaseIdx = FluidSystem::wPhaseIdx,
69 nPhaseIdx = FluidSystem::nPhaseIdx,
70 gPhaseIdx = FluidSystem::gPhaseIdx,
71
72 wCompIdx = FluidSystem::wCompIdx,
73 nCompIdx = FluidSystem::nCompIdx,
74 gCompIdx = FluidSystem::gCompIdx
75 };
76
77public:
78
79 using ParentType::ParentType;
80
92 NumEqVector computeStorage(const Problem& problem,
93 const SubControlVolume& scv,
94 const VolumeVariables& volVars) const
95 {
96 NumEqVector storage(0.0);
97
98 // compute storage term of all components within all phases
99 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
100 {
101 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
102 {
103 auto eqIdx = Indices::conti0EqIdx + compIdx;
104 storage[eqIdx] += volVars.porosity()
105 * volVars.saturation(phaseIdx)
106 * volVars.molarDensity(phaseIdx)
107 * volVars.moleFraction(phaseIdx, compIdx);
108 }
109
110 // The energy storage in the fluid phase with index phaseIdx
111 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
112 }
113
114 // The energy storage in the solid matrix
115 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
116
117 return storage;
118 }
119
131 NumEqVector computeFlux(const Problem& problem,
132 const Element& element,
133 const FVElementGeometry& fvGeometry,
134 const ElementVolumeVariables& elemVolVars,
135 const SubControlVolumeFace& scvf,
136 const ElementFluxVariablesCache& elemFluxVarsCache) const
137 {
138 FluxVariables fluxVars;
139 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
140 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
141
142 // get upwind weights into local scope
143 NumEqVector flux(0.0);
144
145 // advective fluxes
146 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
147 {
148 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
149 {
150 auto upwindTerm = [phaseIdx, compIdx](const VolumeVariables& volVars)
151 { return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
152
153 // get equation index
154 auto eqIdx = Indices::conti0EqIdx + compIdx;
155 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
156 }
157
158 // Add advective phase energy fluxes. For isothermal model the contribution is zero.
159 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
160 }
161
162 // Add diffusive energy fluxes. For isothermal model the contribution is zero.
163 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
164
165 // diffusive fluxes
166 const auto diffusionFluxesWPhase = fluxVars.molecularDiffusionFlux(wPhaseIdx);
167 Scalar jGW = diffusionFluxesWPhase[gCompIdx];
168 Scalar jNW = diffusionFluxesWPhase[nCompIdx];
169 Scalar jWW = -(jGW+jNW);
170
171 //check for the reference system and adapt units of the diffusive flux accordingly.
172 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
173 {
174 jGW /= FluidSystem::molarMass(gCompIdx);
175 jNW /= FluidSystem::molarMass(nCompIdx);
176 jWW /= FluidSystem::molarMass(wCompIdx);
177 }
178
179 const auto diffusionFluxesGPhase = fluxVars.molecularDiffusionFlux(gPhaseIdx);
180 Scalar jWG = diffusionFluxesGPhase[wCompIdx];
181 Scalar jNG = diffusionFluxesGPhase[nCompIdx];
182 Scalar jGG = -(jWG+jNG);
183
184 //check for the reference system and adapt units of the diffusive flux accordingly.
185 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
186 {
187 jWG /= FluidSystem::molarMass(wCompIdx);
188 jNG /= FluidSystem::molarMass(nCompIdx);
189 jGG /= FluidSystem::molarMass(gCompIdx);
190 }
191
192 // At the moment we do not consider diffusion in the NAPL phase
193 const Scalar jWN = 0.0;
194 const Scalar jGN = 0.0;
195 const Scalar jNN = 0.0;
196
197 flux[contiWEqIdx] += jWW+jWG+jWN;
198 flux[contiNEqIdx] += jNW+jNG+jNN;
199 flux[contiGEqIdx] += jGW+jGG+jGN;
200
201 return flux;
202 }
203};
204
205} // end namespace Dumux
206
207#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...
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 Jacobian matrix for problems using the three-phase three-component fu...
Definition: porousmediumflow/3p3c/localresidual.hh:42
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/3p3c/localresidual.hh:92
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/3p3c/localresidual.hh:131
Declares all properties used in Dumux.