Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup ThreePWaterOilModel | ||
10 | * \brief Element-wise calculation of the Jacobian matrix for problems | ||
11 | * using the three-phase three-component fully implicit model. | ||
12 | */ | ||
13 | |||
14 | #ifndef DUMUX_3P2CNI_LOCAL_RESIDUAL_HH | ||
15 | #define DUMUX_3P2CNI_LOCAL_RESIDUAL_HH | ||
16 | |||
17 | #include <dumux/common/properties.hh> | ||
18 | #include <dumux/common/numeqvector.hh> | ||
19 | #include <dumux/flux/referencesystemformulation.hh> | ||
20 | |||
21 | namespace Dumux { | ||
22 | /*! | ||
23 | * \ingroup ThreePWaterOilModel | ||
24 | * \brief Element-wise calculation of the local residual for problems | ||
25 | * using the ThreePWaterOil fully implicit model. | ||
26 | */ | ||
27 | template<class TypeTag> | ||
28 | class ThreePWaterOilLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual> | ||
29 | { | ||
30 | protected: | ||
31 | using ParentType = GetPropType<TypeTag, Properties::BaseLocalResidual>; | ||
32 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
33 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
34 | using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView; | ||
35 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
36 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
37 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
38 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
39 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
40 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
41 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
42 | using Element = typename GridView::template Codim<0>::Entity; | ||
43 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
44 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
45 | using EnergyLocalResidual = GetPropType<TypeTag, Properties::EnergyLocalResidual>; | ||
46 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
47 | |||
48 | enum { | ||
49 | numPhases = GetPropType<TypeTag, Properties::ModelTraits>::numFluidPhases(), | ||
50 | numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents(), | ||
51 | |||
52 | conti0EqIdx = Indices::conti0EqIdx,//!< Index of the mass conservation equation for the water component | ||
53 | conti1EqIdx = conti0EqIdx + 1,//!< Index of the mass conservation equation for the contaminant component | ||
54 | |||
55 | wPhaseIdx = FluidSystem::wPhaseIdx, | ||
56 | nPhaseIdx = FluidSystem::nPhaseIdx, | ||
57 | gPhaseIdx = FluidSystem::gPhaseIdx, | ||
58 | |||
59 | wCompIdx = FluidSystem::wCompIdx, | ||
60 | nCompIdx = FluidSystem::nCompIdx, | ||
61 | }; | ||
62 | |||
63 | //! Property that defines whether mole or mass fractions are used | ||
64 | static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>(); | ||
65 | public: | ||
66 |
2/4✓ Branch 1 taken 153600 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 153600 times.
✗ Branch 5 not taken.
|
307200 | using ParentType::ParentType; |
67 | |||
68 | /*! | ||
69 | * \brief Evaluates the amount of all conservation quantities | ||
70 | * (e.g. phase mass) within a sub-control volume. | ||
71 | * | ||
72 | * The result should be averaged over the volume (e.g. phase mass | ||
73 | * inside a sub control volume divided by the volume) | ||
74 | * | ||
75 | * \param problem The problem | ||
76 | * \param scv The sub-control-volume | ||
77 | * \param volVars The volume variables | ||
78 | */ | ||
79 | 15974400 | NumEqVector computeStorage(const Problem& problem, | |
80 | const SubControlVolume& scv, | ||
81 | const VolumeVariables& volVars) const | ||
82 | { | ||
83 | 15974400 | NumEqVector storage(0.0); | |
84 | |||
85 | ✗ | const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx) | |
86 | 191692800 | { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); }; | |
87 | |||
88 | ✗ | const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx) | |
89 | 191692800 | { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); }; | |
90 | |||
91 | // compute storage term of all components within all phases | ||
92 |
2/2✓ Branch 0 taken 47923200 times.
✓ Branch 1 taken 15974400 times.
|
63897600 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
93 | { | ||
94 |
2/2✓ Branch 0 taken 95846400 times.
✓ Branch 1 taken 47923200 times.
|
143769600 | for (int compIdx = 0; compIdx < numComponents; ++compIdx) |
95 | { | ||
96 | 95846400 | int eqIdx = (compIdx == wCompIdx) ? conti0EqIdx : conti1EqIdx; | |
97 | 95846400 | storage[eqIdx] += volVars.porosity() | |
98 | 191692800 | * volVars.saturation(phaseIdx) | |
99 | 95846400 | * massOrMoleDensity(volVars, phaseIdx) | |
100 | 191692800 | * massOrMoleFraction(volVars, phaseIdx, compIdx); | |
101 | } | ||
102 | |||
103 | //! The energy storage in the fluid phase with index phaseIdx | ||
104 | 95846400 | EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx); | |
105 | } | ||
106 | |||
107 | //! The energy storage in the solid matrix | ||
108 | 31948800 | EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars); | |
109 | |||
110 | 15974400 | return storage; | |
111 | } | ||
112 | |||
113 | /*! | ||
114 | * \brief Evaluates the total flux of all conservation quantities | ||
115 | * over a face of a sub-control volume. | ||
116 | * | ||
117 | * \param problem The problem | ||
118 | * \param element The element | ||
119 | * \param fvGeometry The finite volume element geometry | ||
120 | * \param elemVolVars The element volume variables | ||
121 | * \param scvf The sub control volume face | ||
122 | * \param elemFluxVarsCache The element flux variables cache | ||
123 | */ | ||
124 | 7987200 | NumEqVector computeFlux(const Problem& problem, | |
125 | const Element& element, | ||
126 | const FVElementGeometry& fvGeometry, | ||
127 | const ElementVolumeVariables& elemVolVars, | ||
128 | const SubControlVolumeFace& scvf, | ||
129 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
130 | { | ||
131 | 7987200 | FluxVariables fluxVars; | |
132 | 7987200 | fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache); | |
133 | |||
134 | // get upwind weights into local scope | ||
135 | 7987200 | NumEqVector flux(0.0); | |
136 | |||
137 | ✗ | const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx) | |
138 | 191692800 | { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); }; | |
139 | |||
140 | ✗ | const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx) | |
141 | 191692800 | { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); }; | |
142 | |||
143 | // advective fluxes | ||
144 |
2/2✓ Branch 0 taken 23961600 times.
✓ Branch 1 taken 7987200 times.
|
31948800 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
145 | { | ||
146 |
2/2✓ Branch 0 taken 47923200 times.
✓ Branch 1 taken 23961600 times.
|
71884800 | for (int compIdx = 0; compIdx < numComponents; ++compIdx) |
147 | { | ||
148 | 335462400 | const auto upwindTerm = [&massOrMoleDensity, &massOrMoleFraction, phaseIdx, compIdx] (const auto& volVars) | |
149 | 383385600 | { return massOrMoleDensity(volVars, phaseIdx)*massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); }; | |
150 | |||
151 | // get equation index | ||
152 | 47923200 | auto eqIdx = conti0EqIdx + compIdx; | |
153 | 47923200 | flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm); | |
154 | } | ||
155 | |||
156 | //! Add advective phase energy fluxes. For isothermal model the contribution is zero. | ||
157 | 23961600 | EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx); | |
158 | } | ||
159 | |||
160 | //! Add diffusive energy fluxes. For isothermal model the contribution is zero. | ||
161 | 7987200 | EnergyLocalResidual::heatConductionFlux(flux, fluxVars); | |
162 | |||
163 | // diffusive fluxes | ||
164 | static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation(); | ||
165 | 7987200 | const auto diffusionFluxesWPhase = fluxVars.molecularDiffusionFlux(wPhaseIdx); | |
166 | 7987200 | Scalar jNW = diffusionFluxesWPhase[nCompIdx]; | |
167 | 7987200 | Scalar jWW = -jNW; | |
168 | // check for the reference system and adapt units of the diffusive flux accordingly. | ||
169 | if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged) | ||
170 | { | ||
171 | 7987200 | jNW /= FluidSystem::molarMass(nCompIdx); | |
172 | 7987200 | jWW /= FluidSystem::molarMass(wCompIdx); | |
173 | } | ||
174 | |||
175 | 7987200 | const auto diffusionFluxesGPhase = fluxVars.molecularDiffusionFlux(gPhaseIdx); | |
176 | 7987200 | Scalar jWG = diffusionFluxesGPhase[wCompIdx]; | |
177 | 7987200 | Scalar jNG = -jWG; | |
178 | // check for the reference system and adapt units of the diffusive flux accordingly. | ||
179 | if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged) | ||
180 | { | ||
181 | 7987200 | jWG /= FluidSystem::molarMass(wCompIdx); | |
182 | 7987200 | jNG /= FluidSystem::molarMass(nCompIdx); | |
183 | } | ||
184 | |||
185 | 7987200 | const auto diffusionFluxesNPhase = fluxVars.molecularDiffusionFlux(nPhaseIdx); | |
186 | 7987200 | Scalar jWN = diffusionFluxesNPhase[wCompIdx]; | |
187 | 7987200 | Scalar jNN = -jWN; | |
188 | // check for the reference system and adapt units of the diffusive flux accordingly. | ||
189 | if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged) | ||
190 | { | ||
191 | 7987200 | jWN /= FluidSystem::molarMass(wCompIdx); | |
192 | 7987200 | jNN /= FluidSystem::molarMass(nCompIdx); | |
193 | } | ||
194 | |||
195 | 15974400 | flux[conti0EqIdx] += jWW+jWG+jWN; | |
196 | 15974400 | flux[conti1EqIdx] += jNW+jNG+jNN; | |
197 | |||
198 | 7987200 | return flux; | |
199 | } | ||
200 | }; | ||
201 | |||
202 | } // end namespace Dumux | ||
203 | |||
204 | #endif | ||
205 |