My Project
Loading...
Searching...
No Matches
MultisegmentWellPrimaryVariables.hpp
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
27#include <opm/simulators/wells/MultisegmentWellEquations.hpp>
28#include <opm/input/eclipse/Schedule/SummaryState.hpp>
29
30#include <array>
31#include <cstddef>
32#include <vector>
33
34namespace Opm
35{
36
37class DeferredLogger;
38template<class Scalar> class MultisegmentWellGeneric;
39template<class FluidSystem, class Indices, class Scalar> class WellInterfaceIndices;
40class WellState;
41
42template<class FluidSystem, class Indices, class Scalar>
44{
45public:
46 // TODO: for now, not considering the polymer, solvent and so on to simplify the development process.
47
48 // TODO: we need to have order for the primary variables and also the order for the well equations.
49 // sometimes, they are similar, while sometimes, they can have very different forms.
50
51 // Table showing the primary variable indices, depending on what phases are present:
52 //
53 // WOG OG WG WO W/O/G (single phase)
54 // WQTotal 0 0 0 0 0
55 // WFrac 1 -1000 1 1 -1000
56 // GFrac 2 1 -1000 -1000 -1000
57 // Spres 3 2 2 2 1
58
59 static constexpr bool has_water = (Indices::waterSwitchIdx >= 0);
60 static constexpr bool has_gas = (Indices::compositionSwitchIdx >= 0);
61 static constexpr bool has_oil = (Indices::numPhases - has_gas - has_water) > 0;
62
63 // In the implementation, one should use has_wfrac_variable
64 // rather than has_water to check if you should do something
65 // with the variable at the WFrac location, similar for GFrac.
66 static constexpr bool has_wfrac_variable = has_water && Indices::numPhases > 1;
67 static constexpr bool has_gfrac_variable = has_gas && has_oil;
68
69 static constexpr int WQTotal = 0;
70 static constexpr int WFrac = has_wfrac_variable ? 1 : -1000;
71 static constexpr int GFrac = has_gfrac_variable ? has_wfrac_variable + 1 : -1000;
72 static constexpr int SPres = has_wfrac_variable + has_gfrac_variable + 1;
73
74 // the number of well equations TODO: it should have a more general strategy for it
75 static constexpr int numWellEq = Indices::numPhases + 1;
76
77 using EvalWell = DenseAd::Evaluation<double, /*size=*/Indices::numEq + numWellEq>;
78
80 using BVectorWell = typename Equations::BVectorWell;
81
83 : well_(well)
84 {}
85
87 void resize(const int numSegments);
88
90 void init();
91
93 void update(const WellState& well_state, const bool stop_or_zero_rate_target);
94
96 void updateNewton(const BVectorWell& dwells,
97 const double relaxation_factor,
98 const double DFLimit,
99 const bool stop_or_zero_rate_target,
100 const double max_pressure_change);
101
104 const double rho,
105 const bool stop_or_zero_rate_target,
106 WellState& well_state,
109
112 EvalWell volumeFractionScaled(const int seg,
113 const int compIdx) const;
114
117 EvalWell surfaceVolumeFraction(const int seg,
118 const int compIdx) const;
119
121 EvalWell getSegmentRateUpwinding(const int seg,
122 const int seg_upwind,
123 const std::size_t comp_idx) const;
124
126 EvalWell getBhp() const;
127
129 EvalWell getSegmentPressure(const int seg) const;
130
132 EvalWell getSegmentRate(const int seg,
133 const int comp_idx) const;
134
136 EvalWell getQs(const int comp_idx) const;
137
139 EvalWell getWQTotal() const;
140
142 const std::array<EvalWell, numWellEq>& eval(const int idx) const
143 { return evaluation_[idx]; }
144
146 const std::array<Scalar, numWellEq>& value(const int idx) const
147 { return value_[idx]; }
148
150 void setValue(const int idx, const std::array<Scalar, numWellEq>& val)
151 { value_[idx] = val; }
152
153private:
155 void processFractions(const int seg);
156
158 EvalWell volumeFraction(const int seg,
159 const unsigned compIdx) const;
160
163 std::vector<std::array<double, numWellEq>> value_;
164
167 std::vector<std::array<EvalWell, numWellEq>> evaluation_;
168
170};
171
172}
173
174#endif // OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
Definition AquiferInterface.hpp:35
Definition DeferredLogger.hpp:57
Definition MultisegmentWellPrimaryVariables.hpp:44
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
Definition MultisegmentWellPrimaryVariables.cpp:618
void resize(const int numSegments)
Resize values and evaluations.
Definition MultisegmentWellPrimaryVariables.cpp:47
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:150
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.
Definition MultisegmentWellPrimaryVariables.cpp:156
EvalWell getSegmentPressure(const int seg) const
Get pressure for a segment.
Definition MultisegmentWellPrimaryVariables.cpp:593
EvalWell getBhp() const
Get bottomhole pressure.
Definition MultisegmentWellPrimaryVariables.cpp:601
EvalWell getWQTotal() const
Get WQTotal.
Definition MultisegmentWellPrimaryVariables.cpp:626
EvalWell volumeFractionScaled(const int seg, const int compIdx) const
Returns scaled volume fraction for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:520
const std::array< Scalar, numWellEq > & value(const int idx) const
Returns a value array.
Definition MultisegmentWellPrimaryVariables.hpp:146
void update(const WellState &well_state, const bool stop_or_zero_rate_target)
Copy values from well state.
Definition MultisegmentWellPrimaryVariables.cpp:68
void copyToWellState(const MultisegmentWellGeneric< Scalar > &mswell, const double rho, const bool stop_or_zero_rate_target, WellState &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Copy values to well state.
Definition MultisegmentWellPrimaryVariables.cpp:214
void init()
Initialize evaluations from values.
Definition MultisegmentWellPrimaryVariables.cpp:55
EvalWell surfaceVolumeFraction(const int seg, const int compIdx) const
Returns surface volume fraction for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:537
const std::array< EvalWell, numWellEq > & eval(const int idx) const
Returns a const ref to an array of evaluations.
Definition MultisegmentWellPrimaryVariables.hpp:142
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.
Definition MultisegmentWellPrimaryVariables.cpp:553
EvalWell getSegmentRate(const int seg, const int comp_idx) const
Get rate for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:609
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:60
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27