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