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 ThreePThreeCModel | ||
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_3P3C_LOCAL_RESIDUAL_HH | ||
15 | #define DUMUX_3P3C_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 | /*! | ||
25 | * \ingroup ThreePThreeCModel | ||
26 | * \brief Element-wise calculation of the Jacobian matrix for problems | ||
27 | * using the three-phase three-component fully implicit model. | ||
28 | */ | ||
29 | template<class TypeTag> | ||
30 | class ThreePThreeCLocalResidual | ||
31 | : public DiscretizationDefaultLocalOperator<TypeTag> | ||
32 | { | ||
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 | contiWEqIdx = Indices::conti0EqIdx + FluidSystem::wPhaseIdx,//!< index of the mass conservation equation for the water component | ||
56 | contiNEqIdx = Indices::conti0EqIdx + FluidSystem::nPhaseIdx,//!< index of the mass conservation equation for the contaminant component | ||
57 | contiGEqIdx = Indices::conti0EqIdx + FluidSystem::gPhaseIdx,//!< index of the mass conservation equation for the gas component | ||
58 | |||
59 | wPhaseIdx = FluidSystem::wPhaseIdx, | ||
60 | nPhaseIdx = FluidSystem::nPhaseIdx, | ||
61 | gPhaseIdx = FluidSystem::gPhaseIdx, | ||
62 | |||
63 | wCompIdx = FluidSystem::wCompIdx, | ||
64 | nCompIdx = FluidSystem::nCompIdx, | ||
65 | gCompIdx = FluidSystem::gCompIdx | ||
66 | }; | ||
67 | |||
68 | public: | ||
69 | |||
70 |
1/2✓ Branch 1 taken 134240 times.
✗ Branch 2 not taken.
|
134240 | using ParentType::ParentType; |
71 | |||
72 | /*! | ||
73 | * \brief Evaluates the amount of all conservation quantities | ||
74 | * (e.g. phase mass) within a sub-control volume. | ||
75 | * | ||
76 | * The result should be averaged over the volume (e.g. phase mass | ||
77 | * inside a sub control volume divided by the volume) | ||
78 | * | ||
79 | * \param problem The problem | ||
80 | * \param scv The sub control volume | ||
81 | * \param volVars The volume variables | ||
82 | */ | ||
83 | 23085760 | NumEqVector computeStorage(const Problem& problem, | |
84 | const SubControlVolume& scv, | ||
85 | const VolumeVariables& volVars) const | ||
86 | { | ||
87 | 23085760 | NumEqVector storage(0.0); | |
88 | |||
89 | // compute storage term of all components within all phases | ||
90 |
2/2✓ Branch 0 taken 69257280 times.
✓ Branch 1 taken 23085760 times.
|
92343040 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
91 | { | ||
92 |
2/2✓ Branch 0 taken 207771840 times.
✓ Branch 1 taken 69257280 times.
|
277029120 | for (int compIdx = 0; compIdx < numComponents; ++compIdx) |
93 | { | ||
94 | 207771840 | auto eqIdx = Indices::conti0EqIdx + compIdx; | |
95 | 207771840 | storage[eqIdx] += volVars.porosity() | |
96 | 207771840 | * volVars.saturation(phaseIdx) | |
97 | 207771840 | * volVars.molarDensity(phaseIdx) | |
98 | 207771840 | * volVars.moleFraction(phaseIdx, compIdx); | |
99 | } | ||
100 | |||
101 | // The energy storage in the fluid phase with index phaseIdx | ||
102 | 66709440 | EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx); | |
103 | } | ||
104 | |||
105 | // The energy storage in the solid matrix | ||
106 | 22236480 | EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars); | |
107 | |||
108 | 23085760 | return storage; | |
109 | } | ||
110 | |||
111 | /*! | ||
112 | * \brief Evaluates the total flux of all conservation quantities | ||
113 | * over a face of a sub-control volume. | ||
114 | * | ||
115 | * \param problem The problem | ||
116 | * \param element The element | ||
117 | * \param fvGeometry The finite volume element geometry | ||
118 | * \param elemVolVars The element volume variables | ||
119 | * \param scvf The sub control volume face | ||
120 | * \param elemFluxVarsCache The element flux variables cache | ||
121 | */ | ||
122 | 13474215 | NumEqVector computeFlux(const Problem& problem, | |
123 | const Element& element, | ||
124 | const FVElementGeometry& fvGeometry, | ||
125 | const ElementVolumeVariables& elemVolVars, | ||
126 | const SubControlVolumeFace& scvf, | ||
127 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
128 | { | ||
129 | 13474215 | FluxVariables fluxVars; | |
130 | 13474215 | fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache); | |
131 | static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation(); | ||
132 | |||
133 | // get upwind weights into local scope | ||
134 | 13474215 | NumEqVector flux(0.0); | |
135 | |||
136 | // advective fluxes | ||
137 |
2/2✓ Branch 0 taken 40422645 times.
✓ Branch 1 taken 13474215 times.
|
53896860 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
138 | { | ||
139 |
2/2✓ Branch 0 taken 121267935 times.
✓ Branch 1 taken 40422645 times.
|
161690580 | for (int compIdx = 0; compIdx < numComponents; ++compIdx) |
140 | { | ||
141 | 242535870 | auto upwindTerm = [phaseIdx, compIdx](const VolumeVariables& volVars) | |
142 | 121267935 | { return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); }; | |
143 | |||
144 | // get equation index | ||
145 | 121267935 | auto eqIdx = Indices::conti0EqIdx + compIdx; | |
146 | 121267935 | flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm); | |
147 | } | ||
148 | |||
149 | // Add advective phase energy fluxes. For isothermal model the contribution is zero. | ||
150 | 38708901 | EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx); | |
151 | } | ||
152 | |||
153 | // Add diffusive energy fluxes. For isothermal model the contribution is zero. | ||
154 | 12902967 | EnergyLocalResidual::heatConductionFlux(flux, fluxVars); | |
155 | |||
156 | // diffusive fluxes | ||
157 | 13474215 | const auto diffusionFluxesWPhase = fluxVars.molecularDiffusionFlux(wPhaseIdx); | |
158 | 13474215 | Scalar jGW = diffusionFluxesWPhase[gCompIdx]; | |
159 | 13474215 | Scalar jNW = diffusionFluxesWPhase[nCompIdx]; | |
160 | 13474215 | Scalar jWW = -(jGW+jNW); | |
161 | |||
162 | //check for the reference system and adapt units of the diffusive flux accordingly. | ||
163 | if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged) | ||
164 | { | ||
165 | 13474215 | jGW /= FluidSystem::molarMass(gCompIdx); | |
166 | 13474215 | jNW /= FluidSystem::molarMass(nCompIdx); | |
167 | 13474215 | jWW /= FluidSystem::molarMass(wCompIdx); | |
168 | } | ||
169 | |||
170 | 13474215 | const auto diffusionFluxesGPhase = fluxVars.molecularDiffusionFlux(gPhaseIdx); | |
171 | 13474215 | Scalar jWG = diffusionFluxesGPhase[wCompIdx]; | |
172 | 13474215 | Scalar jNG = diffusionFluxesGPhase[nCompIdx]; | |
173 | 13474215 | Scalar jGG = -(jWG+jNG); | |
174 | |||
175 | //check for the reference system and adapt units of the diffusive flux accordingly. | ||
176 | if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged) | ||
177 | { | ||
178 | 13474215 | jWG /= FluidSystem::molarMass(wCompIdx); | |
179 | 13474215 | jNG /= FluidSystem::molarMass(nCompIdx); | |
180 | 13474215 | jGG /= FluidSystem::molarMass(gCompIdx); | |
181 | } | ||
182 | |||
183 | // At the moment we do not consider diffusion in the NAPL phase | ||
184 | 13474215 | const Scalar jWN = 0.0; | |
185 | 13474215 | const Scalar jGN = 0.0; | |
186 | 13474215 | const Scalar jNN = 0.0; | |
187 | |||
188 | 13474215 | flux[contiWEqIdx] += jWW+jWG+jWN; | |
189 | 13474215 | flux[contiNEqIdx] += jNW+jNG+jNN; | |
190 | 13474215 | flux[contiGEqIdx] += jGW+jGG+jGN; | |
191 | |||
192 | 13474215 | return flux; | |
193 | } | ||
194 | }; | ||
195 | |||
196 | } // end namespace Dumux | ||
197 | |||
198 | #endif | ||
199 |