3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/porousmediumflow/3pwateroil/implicit/problem.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 vesion. *
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_SAGDPROBLEM_HH
26#define DUMUX_SAGDPROBLEM_HH
27
28#include <dune/grid/yaspgrid.hh>
29
31
34
38
39#include "spatialparams.hh"
40
41namespace Dumux {
42
43template <class TypeTag>
44class SagdProblem;
45
46namespace Properties {
47// Create new type tags
48namespace TTag {
49struct Sagd { using InheritsFrom = std::tuple<ThreePWaterOilNI>; };
50struct ThreePWaterOilSagdBox { using InheritsFrom = std::tuple<Sagd, BoxModel>; };
51} // end namespace TTag
52
53// Set the grid type
54template<class TypeTag>
55struct Grid<TypeTag, TTag::Sagd> { using type = Dune::YaspGrid<2>; };
56
57// Set the problem property
58template<class TypeTag>
59struct Problem<TypeTag, TTag::Sagd> { using type = Dumux::SagdProblem<TypeTag>; };
60
61// Set the spatial parameters
62template<class TypeTag>
63struct SpatialParams<TypeTag, TTag::Sagd>
64{
68};
69
70// Set the fluid system
71template<class TypeTag>
72struct FluidSystem<TypeTag, TTag::Sagd>
74
75template<class TypeTag>
76struct OnlyGasPhaseCanDisappear<TypeTag, TTag::Sagd> { static constexpr bool value = true; };
77
78template<class TypeTag>
79struct UseMoles<TypeTag, TTag::Sagd> { static constexpr bool value = true; };
80
81// Set the solid system
82template<class TypeTag>
83struct SolidSystem<TypeTag, TTag::Sagd>
84{
88};
89} // end namespace Properties
90
91
97template <class TypeTag >
98class SagdProblem : public PorousMediumFlowProblem<TypeTag>
99{
106 enum {
107 pressureIdx = Indices::pressureIdx,
108 switch1Idx = Indices::switch1Idx,
109 switch2Idx = Indices::switch2Idx,
110
111 contiWEqIdx = Indices::conti0EqIdx + FluidSystem::wCompIdx,
112 contiNEqIdx = Indices::conti0EqIdx + FluidSystem::nCompIdx,
113 energyEqIdx = Indices::energyEqIdx,
114
115 // phase indices
116 wPhaseIdx = FluidSystem::wPhaseIdx,
117 nPhaseIdx = FluidSystem::nPhaseIdx,
118
119 // phase state
120 wnPhaseOnly = Indices::wnPhaseOnly,
121
122 // world dimension
123 dimWorld = GridView::dimensionworld
124 };
125
129 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
130 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
132 using Element = typename GridView::template Codim<0>::Entity;
133 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
134 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
135 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
136
137public:
138
139 SagdProblem(std::shared_ptr<const GridGeometry> gridGeometry)
140 : ParentType(gridGeometry), pOut_(4e6)
141 {
142 maxDepth_ = 400.0; // [m]
143 FluidSystem::init();
144 totalMassProducedOil_ = 0.0;
145 totalMassProducedWater_ = 0.0;
146 }
147
151 // \{
152
159 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
160 {
161 BoundaryTypes bcTypes;
162 // on bottom
163 if (globalPos[1] < eps_)
164 bcTypes.setAllNeumann();
165
166 // on top
167 else if (globalPos[1] > 40.0 - eps_)
168 bcTypes.setAllNeumann();
169
170 // on bottom other than corners
171 else if (globalPos[0] > 60 - eps_ )
172 bcTypes.setAllDirichlet();
173
174 // on Left
175 else
176 bcTypes.setAllNeumann();
177 return bcTypes;
178 }
179
185 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
186 {
187 return initial_(globalPos);
188 }
189
201 NumEqVector neumann(const Element& element,
202 const FVElementGeometry& fvGeometry,
203 const ElementVolumeVariables& elemVolVars,
204 const ElementFluxVariablesCache& elemFluxVarsCache,
205 const SubControlVolumeFace& scvf) const
206 {
207 NumEqVector values(0.0);
208
209 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
210 const auto& globalPos = insideScv.dofPosition();
211
212 // negative values for injection at injection well
213 if (globalPos[1] > 8.5 - eps_ && globalPos[1] < 9.5 + eps_)
214 {
215 values[contiNEqIdx] = -0.0;
216 values[contiWEqIdx] = -0.193;//*0.5; // (55.5 mol*12.5)/3600 mol/s m = 0.193
217 values[energyEqIdx] = -9132;//*0.5; // J/sec m 9132
218 }
219 else if (globalPos[1] > 2.5 - eps_ && globalPos[1] < 3.5 + eps_) // production well
220 {
221
222 const Scalar elemPressW = elemVolVars[scvf.insideScvIdx()].pressure(wPhaseIdx); //Pressures
223 const Scalar elemPressN = elemVolVars[scvf.insideScvIdx()].pressure(nPhaseIdx);
224
225 const Scalar densityW = elemVolVars[scvf.insideScvIdx()].fluidState().density(wPhaseIdx); //Densities
226 const Scalar densityN = elemVolVars[scvf.insideScvIdx()].fluidState().density(nPhaseIdx);
227
228 const Scalar elemMobW = elemVolVars[scvf.insideScvIdx()].mobility(wPhaseIdx); //Mobilities
229 const Scalar elemMobN = elemVolVars[scvf.insideScvIdx()].mobility(nPhaseIdx);
230
231 const Scalar enthW = elemVolVars[scvf.insideScvIdx()].enthalpy(wPhaseIdx); //Enthalpies
232 const Scalar enthN = elemVolVars[scvf.insideScvIdx()].enthalpy(nPhaseIdx);
233
234 const Scalar wellRadius = 0.50 * 0.3048; // 0.50 ft as specified by SPE9
235
236
237 const Scalar gridHeight_ = 0.5;
238 const Scalar effectiveRadius_ = 0.208 * gridHeight_; //Peaceman's Well Model
239
240 using std::log;
241 // divided by molarMass() of water to convert from kg/m s to mol/m s
242 const Scalar qW = (((2*3.1415*0.5*4e-14)/(log(effectiveRadius_/wellRadius))) *
243 densityW * elemMobW * ( elemPressW-pOut_))/0.01801528;
244 // divided by molarMass() of HeavyOil to convert from kg/m s to mol/m s
245 const Scalar qN = (((2*3.1415*0.5*4e-14)/(log(effectiveRadius_/wellRadius))) *
246 densityN * elemMobN * (elemPressN-pOut_))/0.35;
247
248 Scalar qE;
249 // without cooling:
250 // qE = qW*0.018*enthW + qN*enthN*0.350;
251
252 // with cooling: see Diplomarbeit Stefan Roll, Sept. 2015
253 Scalar wT = elemVolVars[scvf.insideScvIdx()].temperature(); // well temperature
254 if ( wT > 495. )
255 {
256 qE = qW*0.018*enthW + qN*enthN*0.350 + (wT-495.)*5000.; // ~3x injected enthalpy
257 std::cout<< "Cooling now! Extracted enthalpy: " << qE << std::endl;
258 } else {
259 qE = qW*0.018*enthW + qN*enthN*0.350;
260 }
261
262
263 values[contiWEqIdx] = qW;
264 values[contiNEqIdx] = qN;
265 values[energyEqIdx] = qE;
266 massProducedOil_ = qN;
267 massProducedWater_ = qW;
268 }
269 return values;
270 }
271
272 // \}
273
277 // \{
278
284 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
285 {
286 return initial_(globalPos);
287 }
288
289private:
290 // internal method for the initial condition (reused for the
291 // dirichlet conditions!)
292 PrimaryVariables initial_(const GlobalPosition &globalPos) const
293 {
294 PrimaryVariables values(0.0);
295 values.setState(wnPhaseOnly);
296 Scalar densityW = 1000.0;
297 values[pressureIdx] = 101300.0 + (maxDepth_ - globalPos[1])*densityW*9.81;
298
299 values[switch1Idx] = 295.13; // temperature
300 values[switch2Idx] = 0.3; // NAPL saturation
301 return values;
302 }
303
304 Scalar maxDepth_;
305 static constexpr Scalar eps_ = 1e-6;
306 Scalar pIn_;
307 Scalar pOut_;
308 Scalar totalMassProducedOil_;
309 Scalar totalMassProducedWater_;
310
311 mutable Scalar massProducedOil_;
312 mutable Scalar massProducedWater_;
313
314 std::string name_;
315};
316} // end namespace Dumux
317
318#endif
Defines a type tag and some properties for models using the box scheme.
Setting constant fluid properties via the input file.
A compositional fluid system with water and heavy oil components in both the liquid and the gas phase...
A solid phase consisting of a single inert solid component.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:50
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:588
The DUNE grid type.
Definition: common/properties.hh:57
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Property whether to use moles or kg as amount unit for balance equations.
Definition: common/properties.hh:102
The type of the spatial parameters object.
Definition: common/properties.hh:221
The type of the fluid system to use.
Definition: common/properties.hh:223
The type of the solid system to use.
Definition: common/properties.hh:227
reduces the phasestates to threePhases and wnPhaseOnly
Definition: common/properties.hh:264
A component which returns run time specified values for all fluid properties.
Definition: constant.hh:58
A compositional fluid system with water and heavy oil components in both the liquid and the gas phase...
Definition: h2oheavyoil.hh:50
A solid phase consisting of a single inert solid component.
Definition: inertsolidphase.hh:41
Base class for all fully implicit porous media problems.
Definition: dumux/porousmediumflow/problem.hh:39
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:99
NumEqVector neumann(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluates the boundary conditions for a Neumann boundary segment.
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:201
SagdProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:139
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:284
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:185
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Specifies which kind of boundary condition should be used for which equation on a given boundary cont...
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:159
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:49
std::tuple< ThreePWaterOilNI > InheritsFrom
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:49
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:50
std::tuple< Sagd, BoxModel > InheritsFrom
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:50
Dune::YaspGrid< 2 > type
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:55
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:66
GetPropType< TypeTag, Properties::GridGeometry > GridGeometry
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:65
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/3pwateroil/implicit/problem.hh:85
Definition of the spatial parameters for the SAGD problem.
Definition: porousmediumflow/3pwateroil/implicit/spatialparams.hh:47
Adaption of the fully implicit scheme to the three-phase three-component flow model.
Base class for all porous media problems.
Definition of the spatial parameters for the MaxwellStefan problem.