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 ExtendedRichardsModel | ||
10 | * \brief Element-wise calculation of the Jacobian matrix for problems | ||
11 | * using the extended Richards fully implicit models. | ||
12 | */ | ||
13 | |||
14 | #ifndef DUMUX_RICHARDSEXTENDED_LOCAL_RESIDUAL_HH | ||
15 | #define DUMUX_RICHARDSEXTENDED_LOCAL_RESIDUAL_HH | ||
16 | |||
17 | #include <dumux/common/properties.hh> | ||
18 | #include <dumux/common/typetraits/typetraits.hh> | ||
19 | #include <dumux/common/numeqvector.hh> | ||
20 | #include <dumux/discretization/method.hh> | ||
21 | #include <dumux/discretization/extrusion.hh> | ||
22 | #include <dumux/flux/referencesystemformulation.hh> | ||
23 | #include <dumux/porousmediumflow/richards/localresidual.hh> | ||
24 | |||
25 | namespace Dumux { | ||
26 | |||
27 | /*! | ||
28 | * \ingroup ExtendedRichardsModel | ||
29 | * \brief Element-wise calculation of the Jacobian matrix for problems | ||
30 | * using the extended Richards fully implicit models. | ||
31 | */ | ||
32 | template<class TypeTag> | ||
33 | class ExtendedRichardsLocalResidual : public RichardsLocalResidual<TypeTag> | ||
34 | { | ||
35 | using Implementation = GetPropType<TypeTag, Properties::LocalResidual>; | ||
36 | |||
37 | using ParentType = RichardsLocalResidual<TypeTag>; | ||
38 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
39 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
40 | using NumEqVector = Dumux::NumEqVector<GetPropType<TypeTag, Properties::PrimaryVariables>>; | ||
41 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
42 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
43 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
44 | using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; | ||
45 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
46 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
47 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
48 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
49 | using GridView = typename GridGeometry::GridView; | ||
50 | using Element = typename GridView::template Codim<0>::Entity; | ||
51 | using EnergyLocalResidual = GetPropType<TypeTag, Properties::EnergyLocalResidual>; | ||
52 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
53 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
54 | // first index for the mass balance | ||
55 | enum { conti0EqIdx = Indices::conti0EqIdx }; | ||
56 | |||
57 | // phase & component indices | ||
58 | static constexpr auto liquidPhaseIdx = FluidSystem::phase0Idx; | ||
59 | static constexpr auto gasPhaseIdx = FluidSystem::phase1Idx; | ||
60 | static constexpr auto liquidCompIdx = FluidSystem::comp0Idx; | ||
61 | |||
62 | public: | ||
63 |
2/4✓ Branch 1 taken 32700 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 32700 times.
✗ Branch 5 not taken.
|
65400 | using ParentType::ParentType; |
64 | |||
65 | /*! | ||
66 | * \brief Evaluates the rate of change of all conservation | ||
67 | * quantites (e.g. phase mass) within a sub-control | ||
68 | * volume of a finite volume element for the immiscible models. | ||
69 | * \param problem The problem | ||
70 | * \param scv The sub control volume | ||
71 | * \param volVars The current or previous volVars | ||
72 | */ | ||
73 | 1225800 | NumEqVector computeStorage(const Problem& problem, | |
74 | const SubControlVolume& scv, | ||
75 | const VolumeVariables& volVars) const | ||
76 | { | ||
77 | 1225800 | NumEqVector storage = ParentType::computeStorage(problem, scv, volVars); | |
78 | // for extended Richards we consider water in air | ||
79 | 2451600 | storage[conti0EqIdx] += volVars.porosity() | |
80 | 2451600 | * volVars.molarDensity(gasPhaseIdx) | |
81 | 1225800 | * volVars.moleFraction(gasPhaseIdx, liquidCompIdx) | |
82 | 1225800 | * FluidSystem::molarMass(liquidCompIdx) | |
83 | 2451600 | * volVars.saturation(gasPhaseIdx); | |
84 | |||
85 | 1225800 | return storage; | |
86 | } | ||
87 | |||
88 | |||
89 | /*! | ||
90 | * \brief Evaluates the mass flux over a face of a sub control volume. | ||
91 | * | ||
92 | * \param problem The problem | ||
93 | * \param element The current element. | ||
94 | * \param fvGeometry The finite-volume geometry | ||
95 | * \param elemVolVars The volume variables of the current element | ||
96 | * \param scvf The sub control volume face to compute the flux on | ||
97 | * \param elemFluxVarsCache The cache related to flux computation | ||
98 | */ | ||
99 | 925830 | NumEqVector computeFlux(const Problem& problem, | |
100 | const Element& element, | ||
101 | const FVElementGeometry& fvGeometry, | ||
102 | const ElementVolumeVariables& elemVolVars, | ||
103 | const SubControlVolumeFace& scvf, | ||
104 | const ElementFluxVariablesCache& elemFluxVarsCache) const | ||
105 | { | ||
106 | 925830 | FluxVariables fluxVars; | |
107 | 925830 | fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache); | |
108 | 925830 | NumEqVector flux = ParentType::computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache); | |
109 | |||
110 | // for extended Richards we consider water vapor diffusion in air | ||
111 | //check for the reference system and adapt units of the diffusive flux accordingly. | ||
112 | if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged) | ||
113 | 925830 | flux[conti0EqIdx] += fluxVars.molecularDiffusionFlux(gasPhaseIdx)[liquidCompIdx]; | |
114 | else | ||
115 | flux[conti0EqIdx] += fluxVars.molecularDiffusionFlux(gasPhaseIdx)[liquidCompIdx]*FluidSystem::molarMass(liquidCompIdx); | ||
116 | |||
117 | 925830 | return flux; | |
118 | } | ||
119 | }; | ||
120 | |||
121 | } // end namespace Dumux | ||
122 | |||
123 | #endif | ||
124 |