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 TwoPOneCModel | ||
10 | * | ||
11 | * \copydoc Dumux::TwoPOneCLocalResidual | ||
12 | */ | ||
13 | |||
14 | #ifndef DUMUX_2P1C_LOCAL_RESIDUAL_HH | ||
15 | #define DUMUX_2P1C_LOCAL_RESIDUAL_HH | ||
16 | |||
17 | #include <dumux/common/numeqvector.hh> | ||
18 | #include <dumux/porousmediumflow/immiscible/localresidual.hh> | ||
19 | |||
20 | namespace Dumux { | ||
21 | /*! | ||
22 | * \ingroup TwoPOneCModel | ||
23 | * \brief Element-wise calculation of the residual for the fully implicit two-phase one-component flow model. | ||
24 | */ | ||
25 | template<class TypeTag> | ||
26 | class TwoPOneCLocalResidual : public ImmiscibleLocalResidual<TypeTag> | ||
27 | { | ||
28 | using ParentType = ImmiscibleLocalResidual<TypeTag>; | ||
29 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
30 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
31 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
32 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
33 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
34 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
35 | using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView; | ||
36 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
37 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
38 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
39 | using Element = typename GridView::template Codim<0>::Entity; | ||
40 | using EnergyLocalResidual = GetPropType<TypeTag, Properties::EnergyLocalResidual>; | ||
41 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
42 | |||
43 | static const auto numPhases = GetPropType<TypeTag, Properties::ModelTraits>::numFluidPhases(); | ||
44 | |||
45 | public: | ||
46 | //! Use the parent type's constructor | ||
47 |
2/4✓ Branch 1 taken 102600 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 102600 times.
✗ Branch 5 not taken.
|
205200 | using ParentType::ParentType; |
48 | |||
49 | //! Evaluate the storage term within a given scv | ||
50 | 1737600 | NumEqVector computeStorage(const Problem& problem, | |
51 | const SubControlVolume& scv, | ||
52 | const VolumeVariables& volVars) const | ||
53 | { | ||
54 | 1737600 | NumEqVector storage(0.0); | |
55 | // Compute storage term of all components within all phases | ||
56 |
2/2✓ Branch 0 taken 3475200 times.
✓ Branch 1 taken 1737600 times.
|
5212800 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
57 | { | ||
58 | 6950400 | storage[Indices::conti0EqIdx] += | |
59 | volVars.porosity() | ||
60 | 13900800 | * volVars.saturation(phaseIdx) * volVars.density(phaseIdx); | |
61 | |||
62 | 6950400 | EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx); | |
63 | } | ||
64 | |||
65 | // The energy storage in the solid matrix | ||
66 | 3475200 | EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars); | |
67 | |||
68 | 1737600 | return storage; | |
69 | } | ||
70 | |||
71 | //! Evaluate the fluxes over a face of a sub control volume | ||
72 | 2608620 | NumEqVector computeFlux(const Problem& problem, | |
73 | const Element& element, | ||
74 | const FVElementGeometry& fvGeometry, | ||
75 | const ElementVolumeVariables& elemVolVars, | ||
76 | const SubControlVolumeFace& scvf, | ||
77 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
78 | { | ||
79 | 2608620 | FluxVariables fluxVars; | |
80 | 2608620 | fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache); | |
81 | |||
82 | 2608620 | NumEqVector flux; | |
83 |
2/2✓ Branch 0 taken 5217240 times.
✓ Branch 1 taken 2608620 times.
|
7825860 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
84 | { | ||
85 | // The physical quantities for which we perform upwinding | ||
86 | 26086200 | auto upwindTerm = [phaseIdx](const auto& volVars) | |
87 | 31303440 | { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx); }; | |
88 | |||
89 | 5217240 | flux[Indices::conti0EqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm); | |
90 | |||
91 | // Add advective phase energy fluxes. For isothermal model the contribution is zero. | ||
92 | 5217240 | EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx); | |
93 | } | ||
94 | |||
95 | // Add diffusive energy fluxes. For isothermal model the contribution is zero. | ||
96 | 2608620 | EnergyLocalResidual::heatConductionFlux(flux, fluxVars); | |
97 | |||
98 | 2608620 | return flux; | |
99 | } | ||
100 | }; | ||
101 | |||
102 | } // end namespace Dumux | ||
103 | |||
104 | #endif | ||
105 |