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 PorousmediumflowModels | ||
10 | * \brief Element-wise calculation of the residual for problems | ||
11 | * using the n-phase immiscible fully implicit models. | ||
12 | */ | ||
13 | #ifndef DUMUX_IMMISCIBLE_LOCAL_RESIDUAL_HH | ||
14 | #define DUMUX_IMMISCIBLE_LOCAL_RESIDUAL_HH | ||
15 | |||
16 | #include <dumux/common/properties.hh> | ||
17 | #include <dumux/common/numeqvector.hh> | ||
18 | #include <dumux/discretization/defaultlocaloperator.hh> | ||
19 | |||
20 | namespace Dumux { | ||
21 | |||
22 | /*! | ||
23 | * \ingroup PorousmediumflowModels | ||
24 | * \brief Element-wise calculation of the residual for problems | ||
25 | * using the n-phase immiscible fully implicit models. | ||
26 | */ | ||
27 | template<class TypeTag> | ||
28 | class ImmiscibleLocalResidual | ||
29 | : public DiscretizationDefaultLocalOperator<TypeTag> | ||
30 | { | ||
31 | using ThisType = ImmiscibleLocalResidual<TypeTag>; | ||
32 | using ParentType = DiscretizationDefaultLocalOperator<TypeTag>; | ||
33 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
34 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
35 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
36 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
37 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
38 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
39 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
40 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
41 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
42 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
43 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
44 | using GridView = typename GridGeometry::GridView; | ||
45 | using Element = typename GridView::template Codim<0>::Entity; | ||
46 | using EnergyLocalResidual = GetPropType<TypeTag, Properties::EnergyLocalResidual>; | ||
47 | |||
48 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
49 | static constexpr int numPhases = ModelTraits::numFluidPhases(); | ||
50 | static constexpr int conti0EqIdx = ModelTraits::Indices::conti0EqIdx; //!< first index for the mass balance | ||
51 | |||
52 | public: | ||
53 |
5/10✓ Branch 1 taken 11387133 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1048757 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 311728 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 916 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 13489 times.
✗ Branch 14 not taken.
|
14421621 | using ParentType::ParentType; |
54 | |||
55 | /*! | ||
56 | * \brief Evaluate the rate of change of all conservation | ||
57 | * quantites (e.g. phase mass) within a sub-control | ||
58 | * volume of a finite volume element for the immiscible models. | ||
59 | * | ||
60 | * \param problem the problem | ||
61 | * \param scv The sub control volume | ||
62 | * \param volVars The current or previous volVars | ||
63 | */ | ||
64 | 228214016 | NumEqVector computeStorage(const Problem& problem, | |
65 | const SubControlVolume& scv, | ||
66 | const VolumeVariables& volVars) const | ||
67 | { | ||
68 | // partial time derivative of the phase mass | ||
69 | 221719242 | NumEqVector storage; | |
70 |
4/4✓ Branch 0 taken 233653066 times.
✓ Branch 1 taken 115875813 times.
✓ Branch 2 taken 212134858 times.
✓ Branch 3 taken 105843429 times.
|
567890207 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
71 | { | ||
72 | 452282698 | auto eqIdx = conti0EqIdx + phaseIdx; | |
73 | 452282698 | storage[eqIdx] = volVars.porosity() | |
74 | 452282698 | * volVars.density(phaseIdx) | |
75 | 444244884 | * volVars.saturation(phaseIdx); | |
76 | |||
77 | //! The energy storage in the fluid phase with index phaseIdx | ||
78 | 20773504 | EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx); | |
79 | } | ||
80 | |||
81 | //! The energy storage in the solid matrix | ||
82 | 9660032 | EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars); | |
83 | |||
84 | 220823242 | return storage; | |
85 | } | ||
86 | |||
87 | |||
88 | /*! | ||
89 | * \brief Evaluate the mass flux over a face of a sub control volume. | ||
90 | * | ||
91 | * \param problem The problem | ||
92 | * \param element The element | ||
93 | * \param fvGeometry The finite volume geometry context | ||
94 | * \param elemVolVars The volume variables for all flux stencil elements | ||
95 | * \param scvf The sub control volume face to compute the flux on | ||
96 | * \param elemFluxVarsCache The cache related to flux computation | ||
97 | */ | ||
98 | 362872589 | NumEqVector computeFlux(const Problem& problem, | |
99 | const Element& element, | ||
100 | const FVElementGeometry& fvGeometry, | ||
101 | const ElementVolumeVariables& elemVolVars, | ||
102 | const SubControlVolumeFace& scvf, | ||
103 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
104 | { | ||
105 | 362872589 | FluxVariables fluxVars; | |
106 | 362872589 | fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache); | |
107 | |||
108 | 362872589 | NumEqVector flux; | |
109 |
2/2✓ Branch 0 taken 614193173 times.
✓ Branch 1 taken 323986435 times.
|
1016895244 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
110 | { | ||
111 | // the physical quantities for which we perform upwinding | ||
112 | 1309805677 | auto upwindTerm = [phaseIdx](const auto& volVars) | |
113 | 655783022 | { return volVars.density(phaseIdx)*volVars.mobility(phaseIdx); }; | |
114 | |||
115 | 654022655 | auto eqIdx = conti0EqIdx + phaseIdx; | |
116 |
2/4✓ Branch 1 taken 12867680 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12867680 times.
✗ Branch 5 not taken.
|
654022655 | flux[eqIdx] = fluxVars.advectiveFlux(phaseIdx, upwindTerm); |
117 | |||
118 | //! Add advective phase energy fluxes. For isothermal model the contribution is zero. | ||
119 |
1/2✓ Branch 1 taken 12867680 times.
✗ Branch 2 not taken.
|
12891501 | EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx); |
120 | } | ||
121 | |||
122 | //! Add diffusive energy fluxes. For isothermal model the contribution is zero. | ||
123 | 6431445 | EnergyLocalResidual::heatConductionFlux(flux, fluxVars); | |
124 | |||
125 | 362872589 | return flux; | |
126 | } | ||
127 | }; | ||
128 | |||
129 | } // end namespace Dumux | ||
130 | |||
131 | #endif | ||
132 |