MultisegmentWell.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_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_HEADER_INCLUDED
24
27
28namespace Opm {
29
30 class DeferredLogger;
31
32 template<typename TypeTag>
33 class MultisegmentWell : public WellInterface<TypeTag>
34 , public MultisegmentWellEval<GetPropType<TypeTag, Properties::FluidSystem>,
35 GetPropType<TypeTag, Properties::Indices>>
36 {
37 public:
41
42 using typename Base::Simulator;
43 using typename Base::IntensiveQuantities;
44 using typename Base::FluidSystem;
45 using typename Base::ModelParameters;
46 using typename Base::MaterialLaw;
47 using typename Base::Indices;
48 using typename Base::RateConverterType;
49 using typename Base::SparseMatrixAdapter;
50 using typename Base::FluidState;
51
54 using Base::Water;
55 using Base::Oil;
56 using Base::Gas;
57
58 using typename Base::Scalar;
59
61 using typename Base::BVector;
62 using typename Base::Eval;
63
64 using typename MSWEval::Equations;
65 using typename MSWEval::EvalWell;
66 using typename MSWEval::BVectorWell;
67 using MSWEval::SPres;
68 using typename Base::PressureMatrix;
69
70 MultisegmentWell(const Well& well,
71 const ParallelWellInfo<Scalar>& pw_info,
72 const int time_step,
73 const ModelParameters& param,
74 const RateConverterType& rate_converter,
75 const int pvtRegionIdx,
76 const int num_components,
77 const int num_phases,
78 const int index_of_well,
79 const std::vector<PerforationData<Scalar>>& perf_data);
80
81 void init(const PhaseUsage* phase_usage_arg,
82 const std::vector<Scalar>& depth_arg,
83 const Scalar gravity_arg,
84 const int num_cells,
85 const std::vector<Scalar>& B_avg,
86 const bool changed_to_open_this_step) override;
87
88 void initPrimaryVariablesEvaluation() override;
89
91 void updateWellStateWithTarget(const Simulator& simulator,
92 const GroupState<Scalar>& group_state,
93 WellState<Scalar>& well_state,
94 DeferredLogger& deferred_logger) const override;
95
98 const WellState<Scalar>& well_state,
99 const std::vector<Scalar>& B_avg,
100 DeferredLogger& deferred_logger,
101 const bool relax_tolerance) const override;
102
104 void apply(const BVector& x, BVector& Ax) const override;
106 void apply(BVector& r) const override;
107
111 const BVector& x,
112 WellState<Scalar>& well_state,
113 DeferredLogger& deferred_logger) override;
114
116 void computeWellPotentials(const Simulator& simulator,
117 const WellState<Scalar>& well_state,
118 std::vector<Scalar>& well_potentials,
119 DeferredLogger& deferred_logger) override;
120
121 void updatePrimaryVariables(const Simulator& simulator,
122 const WellState<Scalar>& well_state,
123 DeferredLogger& deferred_logger) override;
124
125 void solveEqAndUpdateWellState(const Simulator& simulator,
126 WellState<Scalar>& well_state,
127 DeferredLogger& deferred_logger) override; // const?
128
129 void calculateExplicitQuantities(const Simulator& simulator,
130 const WellState<Scalar>& well_state,
131 DeferredLogger& deferred_logger) override; // should be const?
132
133 void updateIPRImplicit(const Simulator& simulator,
134 WellState<Scalar>& well_state,
135 DeferredLogger& deferred_logger) override;
136
137 void updateProductivityIndex(const Simulator& simulator,
138 const WellProdIndexCalculator<Scalar>& wellPICalc,
139 WellState<Scalar>& well_state,
140 DeferredLogger& deferred_logger) const override;
141
142 Scalar connectionDensity(const int globalConnIdx,
143 const int openConnIdx) const override;
144
145 void addWellContributions(SparseMatrixAdapter& jacobian) const override;
146
148 const BVector& x,
149 const int pressureVarIndex,
150 const bool use_well_weights,
151 const WellState<Scalar>& well_state) const override;
152
153 std::vector<Scalar>
154 computeCurrentWellRates(const Simulator& simulator,
155 DeferredLogger& deferred_logger) const override;
156
157 std::optional<Scalar>
159 const SummaryState& summary_state,
160 const Scalar alq_value,
161 DeferredLogger& deferred_logger) const override;
162
163 std::vector<Scalar> getPrimaryVars() const override;
164
165 int setPrimaryVars(typename std::vector<Scalar>::const_iterator it) override;
166
167 protected:
168 // regularize msw equation
170
171 // the intial amount of fluids in each segment under surface condition
172 std::vector<std::vector<Scalar> > segment_fluid_initial_;
173
174 mutable int debug_cost_counter_ = 0;
175
176 // updating the well_state based on well solution dwells
177 void updateWellState(const Simulator& simulator,
178 const BVectorWell& dwells,
179 WellState<Scalar>& well_state,
180 DeferredLogger& deferred_logger,
181 const Scalar relaxation_factor = 1.0);
182
183 // computing the accumulation term for later use in well mass equations
184 void computeInitialSegmentFluids(const Simulator& simulator);
185
186 // compute the pressure difference between the perforation and cell center
187 void computePerfCellPressDiffs(const Simulator& simulator);
188
189 template<class Value>
190 void computePerfRate(const IntensiveQuantities& int_quants,
191 const std::vector<Value>& mob_perfcells,
192 const std::vector<Scalar>& Tw,
193 const int seg,
194 const int perf,
195 const Value& segment_pressure,
196 const bool& allow_cf,
197 std::vector<Value>& cq_s,
198 Value& perf_press,
199 PerforationRates<Scalar>& perf_rates,
200 DeferredLogger& deferred_logger) const;
201
202 template<class Value>
203 void computePerfRate(const Value& pressure_cell,
204 const Value& rs,
205 const Value& rv,
206 const std::vector<Value>& b_perfcells,
207 const std::vector<Value>& mob_perfcells,
208 const std::vector<Scalar>& Tw,
209 const int perf,
210 const Value& segment_pressure,
211 const Value& segment_density,
212 const bool& allow_cf,
213 const std::vector<Value>& cmix_s,
214 std::vector<Value>& cq_s,
215 Value& perf_press,
216 PerforationRates<Scalar>& perf_rates,
217 DeferredLogger& deferred_logger) const;
218
219 // compute the fluid properties, such as densities, viscosities, and so on, in the segments
220 // They will be treated implicitly, so they need to be of Evaluation type
221 void computeSegmentFluidProperties(const Simulator& simulator,
222 DeferredLogger& deferred_logger);
223
224 // get the mobility for specific perforation
225 template<class Value>
226 void getMobility(const Simulator& simulator,
227 const int perf,
228 std::vector<Value>& mob,
229 DeferredLogger& deferred_logger) const;
230
231 void computeWellRatesAtBhpLimit(const Simulator& simulator,
232 std::vector<Scalar>& well_flux,
233 DeferredLogger& deferred_logger) const;
234
235 void computeWellRatesWithBhp(const Simulator& simulator,
236 const Scalar& bhp,
237 std::vector<Scalar>& well_flux,
238 DeferredLogger& deferred_logger) const override;
239
240 void computeWellRatesWithBhpIterations(const Simulator& simulator,
241 const Scalar& bhp,
242 std::vector<Scalar>& well_flux,
243 DeferredLogger& deferred_logger) const override;
244
245 std::vector<Scalar>
247 const Simulator& simulator,
248 DeferredLogger& deferred_logger) const;
249
250 bool computeWellPotentialsImplicit(const Simulator& simulator,
251 std::vector<Scalar>& well_potentials,
252 DeferredLogger& deferred_logger) const;
253
254 Scalar getRefDensity() const override;
255
256 bool iterateWellEqWithControl(const Simulator& simulator,
257 const double dt,
258 const Well::InjectionControls& inj_controls,
259 const Well::ProductionControls& prod_controls,
260 WellState<Scalar>& well_state,
261 const GroupState<Scalar>& group_state,
262 DeferredLogger& deferred_logger) override;
263
264 bool iterateWellEqWithSwitching(const Simulator& simulator,
265 const double dt,
266 const Well::InjectionControls& inj_controls,
267 const Well::ProductionControls& prod_controls,
268 WellState<Scalar>& well_state,
269 const GroupState<Scalar>& group_state,
270 DeferredLogger& deferred_logger,
271 const bool fixed_control = false,
272 const bool fixed_status = false) override;
273
274 void assembleWellEqWithoutIteration(const Simulator& simulator,
275 const double dt,
276 const Well::InjectionControls& inj_controls,
277 const Well::ProductionControls& prod_controls,
278 WellState<Scalar>& well_state,
279 const GroupState<Scalar>& group_state,
280 DeferredLogger& deferred_logger) override;
281
282 void updateWaterThroughput(const double dt, WellState<Scalar>& well_state) const override;
283
284 EvalWell getSegmentSurfaceVolume(const Simulator& simulator, const int seg_idx) const;
285
286 // turn on crossflow to avoid singular well equations
287 // when the well is banned from cross-flow and the BHP is not properly initialized,
288 // we turn on crossflow to avoid singular well equations. It can result in wrong-signed
289 // well rates, it can cause problem for THP calculation
290 // TODO: looking for better alternative to avoid wrong-signed well rates
291 bool openCrossFlowAvoidSingularity(const Simulator& simulator) const;
292
293 // for a well, when all drawdown are in the wrong direction, then this well will not
294 // be able to produce/inject .
295 bool allDrawDownWrongDirection(const Simulator& simulator) const;
296
297 std::optional<Scalar>
299 const Simulator& ebos_simulator,
300 const SummaryState& summary_state,
301 DeferredLogger& deferred_logger) const;
302
303 std::optional<Scalar>
304 computeBhpAtThpLimitInj(const Simulator& ebos_simulator,
305 const SummaryState& summary_state,
306 DeferredLogger& deferred_logger) const;
307
308 Scalar maxPerfPress(const Simulator& simulator) const;
309
310 // check whether the well is operable under BHP limit with current reservoir condition
312 const Simulator& ebos_simulator,
313 DeferredLogger& deferred_logger) override;
314
315 // check whether the well is operable under THP limit with current reservoir condition
316 void checkOperabilityUnderTHPLimit(const Simulator& ebos_simulator,
317 const WellState<Scalar>& well_state,
318 DeferredLogger& deferred_logger) override;
319
320 // updating the inflow based on the current reservoir condition
321 void updateIPR(const Simulator& ebos_simulator,
322 DeferredLogger& deferred_logger) const override;
323 };
324
325}
326
327#ifndef OPM_MULTISEGMENTWELL_IMPL_HEADER_INCLUDED
329#endif
330
331#endif // OPM_MULTISEGMENTWELL_HEADER_INCLUDED
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupState.hpp:38
Definition: MultisegmentWellEval.hpp:46
MultisegmentWellEquations< Scalar, numWellEq, Indices::numEq > Equations
Definition: MultisegmentWellEval.hpp:54
Definition: MultisegmentWell.hpp:36
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2053
std::vector< Scalar > computeWellPotentialWithTHP(const WellState< Scalar > &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:461
void computeWellPotentials(const Simulator &simulator, const WellState< Scalar > &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: MultisegmentWell_impl.hpp:284
void calculateExplicitQuantities(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:724
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:813
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:835
std::optional< Scalar > computeBhpAtThpLimitProd(const WellState< Scalar > &well_state, const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2036
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1156
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false) override
Definition: MultisegmentWell_impl.hpp:1596
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2216
void solveEqAndUpdateWellState(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:581
ConvergenceReport getWellConvergence(const Simulator &simulator, const WellState< Scalar > &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: MultisegmentWell_impl.hpp:202
void updatePrimaryVariables(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:165
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:741
void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:389
std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2181
void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: MultisegmentWell_impl.hpp:228
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1116
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellState< Scalar > &well_state, DeferredLogger &deferred_logger, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:684
void computePerfRate(const IntensiveQuantities &int_quants, const std::vector< Value > &mob_perfcells, const std::vector< Scalar > &Tw, const int seg, const int perf, const Value &segment_pressure, const bool &allow_cf, std::vector< Value > &cq_s, Value &perf_press, PerforationRates< Scalar > &perf_rates, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1005
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1445
void checkOperabilityUnderBHPLimit(const WellState< Scalar > &well_state, const Simulator &ebos_simulator, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1164
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:1945
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1078
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:263
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: MultisegmentWell_impl.hpp:2234
void computeWellRatesWithBhp(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:345
EvalWell getSegmentSurfaceVolume(const Simulator &simulator, const int seg_idx) const
Definition: MultisegmentWell_impl.hpp:2009
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:1954
void updateWellStateWithTarget(const Simulator &simulator, const GroupState< Scalar > &group_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
updating the well state based the current control mode
Definition: MultisegmentWell_impl.hpp:181
int debug_cost_counter_
Definition: MultisegmentWell.hpp:174
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellState< Scalar > &well_state) const override
Definition: MultisegmentWell_impl.hpp:844
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1789
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1393
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:611
void updateIPR(const Simulator &ebos_simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:1230
void updateWaterThroughput(const double dt, WellState< Scalar > &well_state) const override
Definition: MultisegmentWell_impl.hpp:1998
void initPrimaryVariablesEvaluation() override
Definition: MultisegmentWell_impl.hpp:153
bool computeWellPotentialsImplicit(const Simulator &simulator, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:512
void computeWellRatesAtBhpLimit(const Simulator &simulator, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:329
bool regularize_
Definition: MultisegmentWell.hpp:169
MultisegmentWell(const Well &well, const ParallelWellInfo< Scalar > &pw_info, const int time_step, const ModelParameters &param, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_components, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: MultisegmentWell_impl.hpp:59
Scalar maxPerfPress(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2158
std::vector< std::vector< Scalar > > segment_fluid_initial_
Definition: MultisegmentWell.hpp:172
void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1324
void init(const PhaseUsage *phase_usage_arg, const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const int num_cells, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step) override
Definition: MultisegmentWell_impl.hpp:118
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2104
void computeInitialSegmentFluids(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:666
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:64
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:63
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:65
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:119
Definition: WellInterface.hpp:73
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:78
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:100
BlackOilFluidState< Eval, FluidSystem, has_temperature, has_energy, Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, Indices::numPhases > FluidState
Definition: WellInterface.hpp:132
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:97
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:82
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:83
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
static constexpr bool has_polymer
Definition: WellInterface.hpp:110
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:80
static constexpr bool has_solvent
Definition: WellInterface.hpp:108
BlackoilModelParameters< Scalar > ModelParameters
Definition: WellInterface.hpp:106
typename Base::Eval Eval
Definition: WellInterface.hpp:95
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:81
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:84
Definition: WellProdIndexCalculator.hpp:37
Definition: WellState.hpp:62
Definition: blackoilboundaryratevector.hh:37
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:235
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:147
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: PerforationData.hpp:41
Definition: BlackoilPhases.hpp:46