version 3.10-dev
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// 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_3P3C_LOCAL_RESIDUAL_HH
15#define DUMUX_3P3C_LOCAL_RESIDUAL_HH
16
20
21namespace Dumux
22{
28template<class TypeTag>
29class ThreePThreeCLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
30{
34 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
35 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
36 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
39 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
42 using Element = typename GridView::template Codim<0>::Entity;
43 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
47
48 enum {
51
52 contiWEqIdx = Indices::conti0EqIdx + FluidSystem::wPhaseIdx,
53 contiNEqIdx = Indices::conti0EqIdx + FluidSystem::nPhaseIdx,
54 contiGEqIdx = Indices::conti0EqIdx + FluidSystem::gPhaseIdx,
55
56 wPhaseIdx = FluidSystem::wPhaseIdx,
57 nPhaseIdx = FluidSystem::nPhaseIdx,
58 gPhaseIdx = FluidSystem::gPhaseIdx,
59
60 wCompIdx = FluidSystem::wCompIdx,
61 nCompIdx = FluidSystem::nCompIdx,
62 gCompIdx = FluidSystem::gCompIdx
63 };
64
65public:
66
67 using ParentType::ParentType;
68
80 NumEqVector computeStorage(const Problem& problem,
81 const SubControlVolume& scv,
82 const VolumeVariables& volVars) const
83 {
84 NumEqVector storage(0.0);
85
86 // compute storage term of all components within all phases
87 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
88 {
89 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
90 {
91 auto eqIdx = Indices::conti0EqIdx + compIdx;
92 storage[eqIdx] += volVars.porosity()
93 * volVars.saturation(phaseIdx)
94 * volVars.molarDensity(phaseIdx)
95 * volVars.moleFraction(phaseIdx, compIdx);
96 }
97
98 // The energy storage in the fluid phase with index phaseIdx
99 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx);
100 }
101
102 // The energy storage in the solid matrix
103 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
104
105 return storage;
106 }
107
119 NumEqVector computeFlux(const Problem& problem,
120 const Element& element,
121 const FVElementGeometry& fvGeometry,
122 const ElementVolumeVariables& elemVolVars,
123 const SubControlVolumeFace& scvf,
124 const ElementFluxVariablesCache& elemFluxVarsCache) const
125 {
126 FluxVariables fluxVars;
127 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
128 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
129
130 // get upwind weights into local scope
131 NumEqVector flux(0.0);
132
133 // advective fluxes
134 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
135 {
136 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
137 {
138 auto upwindTerm = [phaseIdx, compIdx](const VolumeVariables& volVars)
139 { return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
140
141 // get equation index
142 auto eqIdx = Indices::conti0EqIdx + compIdx;
143 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
144 }
145
146 // Add advective phase energy fluxes. For isothermal model the contribution is zero.
147 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
148 }
149
150 // Add diffusive energy fluxes. For isothermal model the contribution is zero.
151 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
152
153 // diffusive fluxes
154 const auto diffusionFluxesWPhase = fluxVars.molecularDiffusionFlux(wPhaseIdx);
155 Scalar jGW = diffusionFluxesWPhase[gCompIdx];
156 Scalar jNW = diffusionFluxesWPhase[nCompIdx];
157 Scalar jWW = -(jGW+jNW);
158
159 //check for the reference system and adapt units of the diffusive flux accordingly.
160 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
161 {
162 jGW /= FluidSystem::molarMass(gCompIdx);
163 jNW /= FluidSystem::molarMass(nCompIdx);
164 jWW /= FluidSystem::molarMass(wCompIdx);
165 }
166
167 const auto diffusionFluxesGPhase = fluxVars.molecularDiffusionFlux(gPhaseIdx);
168 Scalar jWG = diffusionFluxesGPhase[wCompIdx];
169 Scalar jNG = diffusionFluxesGPhase[nCompIdx];
170 Scalar jGG = -(jWG+jNG);
171
172 //check for the reference system and adapt units of the diffusive flux accordingly.
173 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
174 {
175 jWG /= FluidSystem::molarMass(wCompIdx);
176 jNG /= FluidSystem::molarMass(nCompIdx);
177 jGG /= FluidSystem::molarMass(gCompIdx);
178 }
179
180 // At the moment we do not consider diffusion in the NAPL phase
181 const Scalar jWN = 0.0;
182 const Scalar jGN = 0.0;
183 const Scalar jNN = 0.0;
184
185 flux[contiWEqIdx] += jWW+jWG+jWN;
186 flux[contiNEqIdx] += jNW+jNG+jNN;
187 flux[contiGEqIdx] += jGW+jGG+jGN;
188
189 return flux;
190 }
191};
192
193} // end namespace Dumux
194
195#endif
Element-wise calculation of the Jacobian matrix for problems using the three-phase three-component fu...
Definition: porousmediumflow/3p3c/localresidual.hh:30
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:80
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:119
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
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...