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 | * \copydoc Dumux::TwoPOneCDarcysLaw | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_2P1C_SPURIOUS_FLUX_BLOCKING_DARCYS_LAW_HH | ||
14 | #define DUMUX_2P1C_SPURIOUS_FLUX_BLOCKING_DARCYS_LAW_HH | ||
15 | |||
16 | #include <dumux/common/math.hh> | ||
17 | #include <dumux/common/parameters.hh> | ||
18 | #include <dumux/common/properties.hh> | ||
19 | #include <dumux/discretization/method.hh> | ||
20 | #include <dumux/flux/darcyslaw.hh> | ||
21 | |||
22 | namespace Dumux { | ||
23 | |||
24 | namespace Properties { | ||
25 | template<class TypeTag, class MyTypeTag> | ||
26 | struct UseBlockingOfSpuriousFlow { using type = UndefinedProperty; }; //!< Determines whether blocking of spurious flow is used or not. | ||
27 | } | ||
28 | |||
29 | /*! | ||
30 | * \ingroup TwoPOneCModel | ||
31 | * \brief Specialization of Darcy's Law for the two-phase one-component model, including a the possibility <BR> | ||
32 | * to block spurious fluxes of cold water into the steam zone, which can improve the model's convergence | ||
33 | * behavior (Gudbjerg et al., 2005) \cite gudbjerg2004. | ||
34 | */ | ||
35 | template <class TypeTag> | ||
36 | class TwoPOneCDarcysLaw : public DarcysLaw<TypeTag> | ||
37 | { | ||
38 | using ParentType = DarcysLaw<TypeTag>; | ||
39 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
40 | using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView; | ||
41 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
42 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
43 | using GridFluxVariablesCache = GetPropType<TypeTag, Properties::GridFluxVariablesCache>; | ||
44 | using ElemFluxVarCache = typename GridFluxVariablesCache::LocalView; | ||
45 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
46 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
47 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
48 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
49 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
50 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
51 | using Element = typename GridView::template Codim<0>::Entity; | ||
52 | using CoordScalar = typename GridView::ctype; | ||
53 | |||
54 | enum { dim = GridView::dimension}; | ||
55 | enum { dimWorld = GridView::dimensionworld}; | ||
56 | |||
57 | // copy some indices for convenience | ||
58 | enum { | ||
59 | // phase indices | ||
60 | liquidPhaseIdx = FluidSystem::liquidPhaseIdx, | ||
61 | gasPhaseIdx = FluidSystem::gasPhaseIdx, | ||
62 | }; | ||
63 | |||
64 | using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>; | ||
65 | using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>; | ||
66 | |||
67 | public: | ||
68 | |||
69 | //! Compute the Darcy flux and use a blocking factor, if specified. | ||
70 | 5217240 | static Scalar flux(const Problem& problem, | |
71 | const Element& element, | ||
72 | const FVElementGeometry& fvGeometry, | ||
73 | const ElementVolumeVariables& elemVolVars, | ||
74 | const SubControlVolumeFace& scvf, | ||
75 | const int phaseIdx, | ||
76 | const ElemFluxVarCache& elemFluxVarCache) | ||
77 | { | ||
78 | 5217240 | const Scalar flux = ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache); | |
79 | |||
80 | // only block wetting-phase (i.e. liquid water) fluxes | ||
81 |
2/2✓ Branch 0 taken 2608620 times.
✓ Branch 1 taken 2608620 times.
|
5217240 | if((!getPropValue<TypeTag, Properties::UseBlockingOfSpuriousFlow>()) || phaseIdx != liquidPhaseIdx) |
82 | return flux; | ||
83 | |||
84 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 612000 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 612000 times.
|
5217240 | const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()]; |
85 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 612000 times.
✓ Branch 2 taken 100687 times.
✓ Branch 3 taken 511313 times.
|
4605240 | const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()]; |
86 | |||
87 |
4/4✓ Branch 0 taken 1083458 times.
✓ Branch 1 taken 1525162 times.
✓ Branch 2 taken 1083458 times.
✓ Branch 3 taken 1525162 times.
|
5217240 | const Scalar factor = std::signbit(flux) ? factor_(outsideVolVars, insideVolVars, phaseIdx) : |
88 | 1525162 | factor_(insideVolVars, outsideVolVars, phaseIdx); | |
89 | |||
90 | 2608620 | return flux * factor; | |
91 | } | ||
92 | |||
93 | private: | ||
94 | /*! | ||
95 | * \brief Calculates the blocking factor which prevents spurious cold water fluxes into the steam zone (Gudbjerg et al., 2005) \cite gudbjerg2004 <BR> | ||
96 | * | ||
97 | * \param up The upstream volume variables | ||
98 | * \param dn The downstream volume variables | ||
99 | * \param phaseIdx The index of the fluid phase | ||
100 | */ | ||
101 | static Scalar factor_(const VolumeVariables &up, const VolumeVariables &dn,const int phaseIdx) | ||
102 | { | ||
103 | 2608620 | Scalar factor = 1.0; | |
104 | |||
105 |
4/4✓ Branch 0 taken 1080554 times.
✓ Branch 1 taken 2904 times.
✓ Branch 2 taken 1520935 times.
✓ Branch 3 taken 4227 times.
|
2608620 | const Scalar tDn = dn.temperature(); //temperature of the downstream SCV (where the cold water is potentially intruding into a steam zone) |
106 |
4/4✓ Branch 0 taken 1080554 times.
✓ Branch 1 taken 2904 times.
✓ Branch 2 taken 1520935 times.
✓ Branch 3 taken 4227 times.
|
2608620 | const Scalar tUp = up.temperature(); //temperature of the upstream SCV |
107 | |||
108 |
4/4✓ Branch 0 taken 1080554 times.
✓ Branch 1 taken 2904 times.
✓ Branch 2 taken 1520935 times.
✓ Branch 3 taken 4227 times.
|
2608620 | const Scalar sgDn = dn.saturation(gasPhaseIdx); //gas phase saturation of the downstream SCV |
109 |
4/4✓ Branch 0 taken 1080554 times.
✓ Branch 1 taken 2904 times.
✓ Branch 2 taken 1520935 times.
✓ Branch 3 taken 4227 times.
|
2608620 | const Scalar sgUp = up.saturation(gasPhaseIdx); //gas phase saturation of the upstream SCV |
110 | |||
111 | 2608620 | bool upIsNotSteam = false; | |
112 | 2608620 | bool downIsSteam = false; | |
113 | 2608620 | bool spuriousFlow = false; | |
114 | |||
115 |
4/4✓ Branch 0 taken 1080554 times.
✓ Branch 1 taken 2904 times.
✓ Branch 2 taken 1520935 times.
✓ Branch 3 taken 4227 times.
|
2608620 | if(sgUp <= 1e-5) |
116 | 2601489 | upIsNotSteam = true; | |
117 | |||
118 |
4/4✓ Branch 0 taken 1980 times.
✓ Branch 1 taken 1081478 times.
✓ Branch 2 taken 1773 times.
✓ Branch 3 taken 1523389 times.
|
2608620 | if(sgDn > 1e-5) |
119 | 3753 | downIsSteam = true; | |
120 | |||
121 |
8/8✓ Branch 0 taken 966 times.
✓ Branch 1 taken 1082492 times.
✓ Branch 2 taken 900 times.
✓ Branch 3 taken 66 times.
✓ Branch 4 taken 759 times.
✓ Branch 5 taken 1524403 times.
✓ Branch 6 taken 693 times.
✓ Branch 7 taken 66 times.
|
2608620 | if(upIsNotSteam && downIsSteam && tDn > tUp && phaseIdx == liquidPhaseIdx) |
122 | 1593 | spuriousFlow = true; | |
123 | |||
124 | if(spuriousFlow) | ||
125 | { | ||
126 | 1593 | Scalar deltaT = tDn - tUp; | |
127 | |||
128 |
4/4✓ Branch 0 taken 150 times.
✓ Branch 1 taken 750 times.
✓ Branch 2 taken 285 times.
✓ Branch 3 taken 408 times.
|
1593 | if((deltaT) > 15 ) |
129 | factor = 0.0 ; | ||
130 | |||
131 | else | ||
132 | 435 | factor = 1-(deltaT/15); | |
133 | } | ||
134 | return factor; | ||
135 | } | ||
136 | }; | ||
137 | |||
138 | } // end namespace Dumux | ||
139 | |||
140 | #endif | ||
141 |