MultisegmentWellPrimaryVariables.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21
22#ifndef OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
24
25#include <opm/material/densead/Evaluation.hpp>
26
28#include <opm/input/eclipse/Schedule/SummaryState.hpp>
29#include <opm/input/eclipse/Units/Units.hpp>
30
31#include <array>
32#include <cstddef>
33#include <vector>
34
35namespace Opm
36{
37
38class DeferredLogger;
39template<class Scalar> class MultisegmentWellGeneric;
40template<class FluidSystem, class Indices> class WellInterfaceIndices;
41template<class Scalar> class WellState;
42
43template<class FluidSystem, class Indices>
45{
46public:
47 // TODO: for now, not considering the polymer, solvent and so on to simplify the development process.
48
49 // TODO: we need to have order for the primary variables and also the order for the well equations.
50 // sometimes, they are similar, while sometimes, they can have very different forms.
51
52 // Table showing the primary variable indices, depending on what phases are present:
53 //
54 // WOG OG WG WO W/O/G (single phase)
55 // WQTotal 0 0 0 0 0
56 // WFrac 1 -1000 1 1 -1000
57 // GFrac 2 1 -1000 -1000 -1000
58 // Spres 3 2 2 2 1
59
60 static constexpr bool has_water = (Indices::waterSwitchIdx >= 0);
61 static constexpr bool has_gas = (Indices::compositionSwitchIdx >= 0);
62 static constexpr bool has_oil = (Indices::numPhases - has_gas - has_water) > 0;
63
64 // In the implementation, one should use has_wfrac_variable
65 // rather than has_water to check if you should do something
66 // with the variable at the WFrac location, similar for GFrac.
67 static constexpr bool has_wfrac_variable = has_water && Indices::numPhases > 1;
68 static constexpr bool has_gfrac_variable = has_gas && has_oil;
69
70 static constexpr int WQTotal = 0;
71 static constexpr int WFrac = has_wfrac_variable ? 1 : -1000;
72 static constexpr int GFrac = has_gfrac_variable ? has_wfrac_variable + 1 : -1000;
73 static constexpr int SPres = has_wfrac_variable + has_gfrac_variable + 1;
74
75 // the number of well equations TODO: it should have a more general strategy for it
76 static constexpr int numWellEq = Indices::numPhases + 1;
77
78 using Scalar = typename FluidSystem::Scalar;
79 using EvalWell = DenseAd::Evaluation<double, /*size=*/Indices::numEq + numWellEq>;
80
83
85 : well_(well)
86 {}
87
89 void resize(const int numSegments);
90
92 void init();
93
95 void update(const WellState<Scalar>& well_state,
96 const bool stop_or_zero_rate_target);
97
99 void updateNewton(const BVectorWell& dwells,
100 const double relaxation_factor,
101 const double DFLimit,
102 const bool stop_or_zero_rate_target,
103 const double max_pressure_change);
104
107 const double rho,
108 const bool stop_or_zero_rate_target,
109 WellState<Scalar>& well_state,
110 const SummaryState& summary_state,
111 DeferredLogger& deferred_logger) const;
112
116 const int compIdx) const;
117
121 const int compIdx) const;
122
125 const int seg_upwind,
126 const std::size_t comp_idx) const;
127
130
132 EvalWell getSegmentPressure(const int seg) const;
133
136 const int comp_idx) const;
137
139 EvalWell getQs(const int comp_idx) const;
140
143
145 const std::array<EvalWell, numWellEq>& eval(const int idx) const
146 { return evaluation_[idx]; }
147
149 const std::array<Scalar, numWellEq>& value(const int idx) const
150 { return value_[idx]; }
151
153 void setValue(const int idx, const std::array<Scalar, numWellEq>& val)
154 { value_[idx] = val; }
155
157 void outputLowLimitPressureSegments(DeferredLogger& deferred_logger) const;
158
159private:
161 void processFractions(const int seg);
162
164 EvalWell volumeFraction(const int seg,
165 const unsigned compIdx) const;
166
169 std::vector<std::array<double, numWellEq>> value_;
170
173 std::vector<std::array<EvalWell, numWellEq>> evaluation_;
174
176
177 static constexpr double bhp_lower_limit = 1. * unit::barsa - 1. * unit::Pascal;
178 static constexpr double seg_pres_lower_limit = 0.;
179};
180
181}
182
183#endif // OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
Definition: DeferredLogger.hpp:57
Dune::BlockVector< VectorBlockWellType > BVectorWell
Definition: MultisegmentWellEquations.hpp:57
Definition: MultisegmentWellGeneric.hpp:42
Definition: MultisegmentWellPrimaryVariables.hpp:45
static constexpr int WFrac
Definition: MultisegmentWellPrimaryVariables.hpp:71
void resize(const int numSegments)
Resize values and evaluations.
DenseAd::Evaluation< double, Indices::numEq+numWellEq > EvalWell
Definition: MultisegmentWellPrimaryVariables.hpp:79
void copyToWellState(const MultisegmentWellGeneric< Scalar > &mswell, const double rho, const bool stop_or_zero_rate_target, WellState< Scalar > &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Copy values to well state.
EvalWell getBhp() const
Get bottomhole pressure.
void updateNewton(const BVectorWell &dwells, const double relaxation_factor, const double DFLimit, const bool stop_or_zero_rate_target, const double max_pressure_change)
Update values from newton update vector.
static constexpr int numWellEq
Definition: MultisegmentWellPrimaryVariables.hpp:76
static constexpr bool has_gfrac_variable
Definition: MultisegmentWellPrimaryVariables.hpp:68
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
void outputLowLimitPressureSegments(DeferredLogger &deferred_logger) const
output the segments with pressure close to lower pressure limit for debugging purpose
typename FluidSystem::Scalar Scalar
Definition: MultisegmentWellPrimaryVariables.hpp:78
EvalWell getWQTotal() const
Get WQTotal.
EvalWell getSegmentPressure(const int seg) const
Get pressure for a segment.
void init()
Initialize evaluations from values.
const std::array< EvalWell, numWellEq > & eval(const int idx) const
Returns a const ref to an array of evaluations.
Definition: MultisegmentWellPrimaryVariables.hpp:145
static constexpr bool has_water
Definition: MultisegmentWellPrimaryVariables.hpp:60
static constexpr bool has_wfrac_variable
Definition: MultisegmentWellPrimaryVariables.hpp:67
EvalWell volumeFractionScaled(const int seg, const int compIdx) const
Returns scaled volume fraction for a component in a segment.
MultisegmentWellPrimaryVariables(const WellInterfaceIndices< FluidSystem, Indices > &well)
Definition: MultisegmentWellPrimaryVariables.hpp:84
static constexpr bool has_oil
Definition: MultisegmentWellPrimaryVariables.hpp:62
void setValue(const int idx, const std::array< Scalar, numWellEq > &val)
Set a value array. Note that this does not also set the corresponding evaluation.
Definition: MultisegmentWellPrimaryVariables.hpp:153
static constexpr int SPres
Definition: MultisegmentWellPrimaryVariables.hpp:73
static constexpr bool has_gas
Definition: MultisegmentWellPrimaryVariables.hpp:61
static constexpr int GFrac
Definition: MultisegmentWellPrimaryVariables.hpp:72
typename Equations::BVectorWell BVectorWell
Definition: MultisegmentWellPrimaryVariables.hpp:82
void update(const WellState< Scalar > &well_state, const bool stop_or_zero_rate_target)
Copy values from well state.
EvalWell getSegmentRateUpwinding(const int seg, const int seg_upwind, const std::size_t comp_idx) const
Returns upwinding rate for a component in a segment.
static constexpr int WQTotal
Definition: MultisegmentWellPrimaryVariables.hpp:70
EvalWell getSegmentRate(const int seg, const int comp_idx) const
Get rate for a component in a segment.
const std::array< Scalar, numWellEq > & value(const int idx) const
Returns a value array.
Definition: MultisegmentWellPrimaryVariables.hpp:149
EvalWell surfaceVolumeFraction(const int seg, const int compIdx) const
Returns surface volume fraction for a component in a segment.
Definition: WellInterfaceIndices.hpp:35
Definition: WellState.hpp:62
Definition: BlackoilPhases.hpp:27