3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/2p1c/darcyslaw.hh
Go to the documentation of this file.
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_2P1C_SPURIOUS_FLUX_BLOCKING_DARCYS_LAW_HH
26#define DUMUX_2P1C_SPURIOUS_FLUX_BLOCKING_DARCYS_LAW_HH
27
28#include <dumux/common/math.hh>
33
34namespace Dumux {
35
36namespace Properties {
37 template<class TypeTag, class MyTypeTag>
39}
40
47template <class TypeTag>
48class TwoPOneCDarcysLaw : public DarcysLaw<TypeTag>
49{
52 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
53 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
54 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
56 using ElemFluxVarCache = typename GridFluxVariablesCache::LocalView;
57 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
63 using Element = typename GridView::template Codim<0>::Entity;
64 using CoordScalar = typename GridView::ctype;
65
66 enum { dim = GridView::dimension};
67 enum { dimWorld = GridView::dimensionworld};
68
69 // copy some indices for convenience
70 enum {
71 // phase indices
72 liquidPhaseIdx = FluidSystem::liquidPhaseIdx,
73 gasPhaseIdx = FluidSystem::gasPhaseIdx,
74 };
75
76 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
77 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
78
79public:
80
82 static Scalar flux(const Problem& problem,
83 const Element& element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const SubControlVolumeFace& scvf,
87 const int phaseIdx,
88 const ElemFluxVarCache& elemFluxVarCache)
89 {
90 const Scalar flux = ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache);
91
92 // only block wetting-phase (i.e. liquid water) fluxes
93 if((!getPropValue<TypeTag, Properties::UseBlockingOfSpuriousFlow>()) || phaseIdx != liquidPhaseIdx)
94 return flux;
95
96 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
97 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
98
99 const Scalar factor = std::signbit(flux) ? factor_(outsideVolVars, insideVolVars, phaseIdx) :
100 factor_(insideVolVars, outsideVolVars, phaseIdx);
101
102 return flux * factor;
103 }
104
105private:
113 static Scalar factor_(const VolumeVariables &up, const VolumeVariables &dn,const int phaseIdx)
114 {
115 Scalar factor = 1.0;
116
117 const Scalar tDn = dn.temperature(); //temperature of the downstream SCV (where the cold water is potentially intruding into a steam zone)
118 const Scalar tUp = up.temperature(); //temperature of the upstream SCV
119
120 const Scalar sgDn = dn.saturation(gasPhaseIdx); //gas phase saturation of the downstream SCV
121 const Scalar sgUp = up.saturation(gasPhaseIdx); //gas phase saturation of the upstream SCV
122
123 bool upIsNotSteam = false;
124 bool downIsSteam = false;
125 bool spuriousFlow = false;
126
127 if(sgUp <= 1e-5)
128 upIsNotSteam = true;
129
130 if(sgDn > 1e-5)
131 downIsSteam = true;
132
133 if(upIsNotSteam && downIsSteam && tDn > tUp && phaseIdx == liquidPhaseIdx)
134 spuriousFlow = true;
135
136 if(spuriousFlow)
137 {
138 Scalar deltaT = tDn - tUp;
139
140 if((deltaT) > 15 )
141 factor = 0.0 ;
142
143 else
144 factor = 1-(deltaT/15);
145 }
146 return factor;
147 }
148};
149
150} // end namespace Dumux
151
152#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
a tag to mark properties as undefined
Definition: propertysystem.hh:35
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:44
Determines whether blocking of spurious flow is used or not.
Definition: porousmediumflow/2p1c/darcyslaw.hh:38
Specialization of Darcy's Law for the two-phase one-component model, including a the possibility to...
Definition: porousmediumflow/2p1c/darcyslaw.hh:49
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElemFluxVarCache &elemFluxVarCache)
Compute the Darcy flux and use a blocking factor, if specified.
Definition: porousmediumflow/2p1c/darcyslaw.hh:82
Declares all properties used in Dumux.
Advective fluxes according to Darcy's law.