BlackoilWellModelGeneric.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 SINTEF ICT, Applied Mathematics.
3 Copyright 2016 - 2017 Statoil ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 IRIS AS
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21*/
22
23#ifndef OPM_BLACKOILWELLMODEL_GENERIC_HEADER_INCLUDED
24#define OPM_BLACKOILWELLMODEL_GENERIC_HEADER_INCLUDED
25
26#include <opm/output/data/GuideRateValue.hpp>
27
28#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
29#include <opm/input/eclipse/Schedule/Well/PAvg.hpp>
30#include <opm/input/eclipse/Schedule/Well/PAvgCalculator.hpp>
31#include <opm/input/eclipse/Schedule/Well/PAvgCalculatorCollection.hpp>
32#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
33
35
42
43#include <cstddef>
44#include <functional>
45#include <map>
46#include <memory>
47#include <optional>
48#include <string>
49#include <unordered_map>
50#include <unordered_set>
51#include <vector>
52
53namespace Opm {
54 class DeferredLogger;
55 class EclipseState;
56 class GasLiftGroupInfo;
57 class GasLiftSingleWellGeneric;
58 class GasLiftWellState;
59 class Group;
60 class GuideRateConfig;
61 class ParallelWellInfo;
62 class RestartValue;
63 class Schedule;
64 struct SimulatorUpdate;
65 class SummaryConfig;
66 class VFPProperties;
67 class WellInterfaceGeneric;
68 template<class Scalar> class WellState;
69} // namespace Opm
70
71namespace Opm { namespace data {
72 struct GroupData;
73 struct GroupGuideRates;
74 class GroupAndNetworkValues;
75 struct NodeData;
76}} // namespace Opm::data
77
78namespace Opm {
79
81template<class Scalar>
83{
84public:
85 // --------- Types ---------
86 using GLiftOptWells = std::map<std::string, std::unique_ptr<GasLiftSingleWellGeneric>>;
87 using GLiftProdWells = std::map<std::string, const WellInterfaceGeneric*>;
88 using GLiftWellStateMap = std::map<std::string, std::unique_ptr<GasLiftWellState>>;
89
91 const SummaryState& summaryState,
92 const EclipseState& eclState,
93 const PhaseUsage& phase_usage,
95
96 virtual ~BlackoilWellModelGeneric() = default;
97
98 int numLocalWells() const;
99 int numLocalWellsEnd() const;
101 int numPhases() const;
102
104 bool wellsActive() const;
105 bool hasWell(const std::string& wname) const;
106
108 bool networkActive() const;
109
110 // whether there exists any multisegment well open on this process
111 bool anyMSWellOpenLocal() const;
112
113 const Well& getWellEcl(const std::string& well_name) const;
114 std::vector<Well> getLocalWells(const int timeStepIdx) const;
115 const Schedule& schedule() const { return schedule_; }
116 const PhaseUsage& phaseUsage() const { return phase_usage_; }
117 const GroupState<Scalar>& groupState() const { return this->active_wgstate_.group_state; }
118 std::vector<const WellInterfaceGeneric*> genericWells() const
119 { return {well_container_generic_.begin(), well_container_generic_.end()}; }
120
121 /*
122 Immutable version of the currently active wellstate.
123 */
125 {
126 return this->active_wgstate_.well_state;
127 }
128
129 /*
130 Mutable version of the currently active wellstate.
131 */
133 {
134 return this->active_wgstate_.well_state;
135 }
136
137 /*
138 Will return the currently active nupcolWellState; must initialize
139 the internal nupcol wellstate with initNupcolWellState() first.
140 */
142 {
143 return this->nupcol_wgstate_.well_state;
144 }
145 GroupState<Scalar>& groupState() { return this->active_wgstate_.group_state; }
146
147 WellTestState& wellTestState() { return this->active_wgstate_.well_test_state; }
148
149 const WellTestState& wellTestState() const { return this->active_wgstate_.well_test_state; }
150
151 Scalar wellPI(const int well_index) const;
152 Scalar wellPI(const std::string& well_name) const;
153
154 void updateEclWells(const int timeStepIdx,
155 const SimulatorUpdate& sim_update,
156 const SummaryState& st);
157
158 void initFromRestartFile(const RestartValue& restartValues,
159 WellTestState wtestState,
160 const std::size_t numCells,
161 bool handle_ms_well);
162
163 void prepareDeserialize(int report_step,
164 const std::size_t numCells,
165 bool handle_ms_well);
166
167 /*
168 Will assign the internal member last_valid_well_state_ to the
169 current value of the this->active_well_state_. The state stored
170 with storeWellState() can then subsequently be recovered with the
171 resetWellState() method.
172 */
174 {
176 }
177
178 data::GroupAndNetworkValues groupAndNetworkData(const int reportStepIdx) const;
179
181 bool hasTHPConstraints() const;
182
184 void updateNetworkActiveState(const int report_step);
185
189 bool needPreStepNetworkRebalance(const int report_step) const;
190
193 bool forceShutWellByName(const std::string& wellname,
194 const double simulation_time);
195
196 const std::vector<PerforationData>& perfData(const int well_idx) const
197 { return well_perf_data_[well_idx]; }
198
199 const Parallel::Communication& comm() const { return comm_; }
200
201 const EclipseState& eclipseState() const { return eclState_; }
202
203 const SummaryState& summaryState() const { return summaryState_; }
204
205 const GuideRate& guideRate() const { return guideRate_; }
206
207 bool reportStepStarts() const { return report_step_starts_; }
208
209 bool shouldBalanceNetwork(const int reportStepIndex,
210 const int iterationIdx) const;
211
212 void updateClosedWellsThisStep(const std::string& well_name) const
213 {
214 this->closed_this_step_.insert(well_name);
215 }
216 bool wasDynamicallyShutThisTimeStep(const std::string& well_name) const;
217
218 template<class Serializer>
219 void serializeOp(Serializer& serializer)
220 {
221 serializer(initial_step_);
222 serializer(report_step_starts_);
223 serializer(last_run_wellpi_);
224 serializer(local_shut_wells_);
225 serializer(closed_this_step_);
226 serializer(guideRate_);
227 serializer(node_pressures_);
228 serializer(prev_inj_multipliers_);
229 serializer(active_wgstate_);
230 serializer(last_valid_wgstate_);
231 serializer(nupcol_wgstate_);
232 serializer(last_glift_opt_time_);
233 serializer(switched_prod_groups_);
234 serializer(switched_inj_groups_);
235 serializer(closed_offending_wells_);
236 }
237
239 {
240 return this->initial_step_ == rhs.initial_step_ &&
242 this->last_run_wellpi_ == rhs.last_run_wellpi_ &&
245 this->node_pressures_ == rhs.node_pressures_ &&
247 this->active_wgstate_ == rhs.active_wgstate_ &&
249 this->nupcol_wgstate_ == rhs.nupcol_wgstate_ &&
254 }
255
256protected:
257 /*
258 The dynamic state of the well model is maintained with an instance
259 of the WellState class. Currently we have
260 three different wellstate instances:
261
262 1. The currently active wellstate is in the active_well_state_
263 member. That is the state which is mutated by the simulator.
264
265 2. In the case timestep fails to converge and we must go back and
266 try again with a smaller timestep we need to recover the last
267 valid wellstate. This is maintained with the
268 last_valid_well_state_ member and the functions
269 commitWellState() and resetWellState().
270
271 3. For the NUPCOL functionality we should either use the
272 currently active wellstate or a wellstate frozen at max
273 nupcol iterations. This is handled with the member
274 nupcol_well_state_ and the initNupcolWellState() function.
275 */
276
277 /*
278 Will return the last good wellstate. This is typcially used when
279 initializing a new report step where the Schedule object might
280 have introduced new wells. The wellstate returned by
281 prevWellState() must have been stored with the commitWellState()
282 function first.
283 */
285 {
286 return this->last_valid_wgstate_.well_state;
287 }
288
290 {
291 return this->last_valid_wgstate_;
292 }
293
294 /*
295 Will store a copy of the input argument well_state in the
296 last_valid_well_state_ member, that state can then be recovered
297 with a subsequent call to resetWellState().
298 */
300 {
301 this->last_valid_wgstate_ = std::move(wgstate);
302 }
303
304 /*
305 Will update the internal variable active_well_state_ to whatever
306 was stored in the last_valid_well_state_ member. This function
307 works in pair with commitWellState() which should be called first.
308 */
310 {
312 }
313
314 /*
315 Will store the current active wellstate in the nupcol_well_state_
316 member. This can then be subsequently retrieved with accessor
317 nupcolWellState().
318 */
320 {
321 this->nupcol_wgstate_ = this->active_wgstate_;
322 }
323
326 std::vector<std::reference_wrapper<ParallelWellInfo>>
327 createLocalParallelWellInfo(const std::vector<Well>& wells);
328
331
332 bool wasDynamicallyShutThisTimeStep(const int well_index) const;
333
334 Scalar updateNetworkPressures(const int reportStepIdx);
335
336 void updateWsolvent(const Group& group,
337 const int reportStepIdx,
339 void setWsolvent(const Group& group,
340 const int reportStepIdx,
341 Scalar wsolvent);
342 virtual void calcRates(const int fipnum,
343 const int pvtreg,
344 const std::vector<Scalar>& production_rates,
345 std::vector<Scalar>& resv_coeff) = 0;
346 virtual void calcInjRates(const int fipnum,
347 const int pvtreg,
348 std::vector<Scalar>& resv_coeff) = 0;
349
350 void assignShutConnections(data::Wells& wsrpt,
351 const int reportStepIndex) const;
352 void assignGroupControl(const Group& group,
353 data::GroupData& gdata) const;
354 void assignGroupValues(const int reportStepIdx,
355 std::map<std::string, data::GroupData>& gvalues) const;
356 void assignNodeValues(std::map<std::string, data::NodeData>& nodevalues,
357 const int reportStepIdx) const;
358
359 void calculateEfficiencyFactors(const int reportStepIdx);
360
361 void checkGconsaleLimits(const Group& group,
362 WellState<Scalar>& well_state,
363 const int reportStepIdx,
364 DeferredLogger& deferred_logger);
365
366 void checkGEconLimits(const Group& group,
367 const double simulation_time,
368 const int report_step_idx,
369 DeferredLogger& deferred_logger);
370
371 bool checkGroupHigherConstraints(const Group& group,
372 DeferredLogger& deferred_logger,
373 const int reportStepIdx);
374
375 void updateAndCommunicateGroupData(const int reportStepIdx,
376 const int iterationIdx);
377
379
381
382 void gliftDebug(const std::string& msg,
383 DeferredLogger& deferred_logger) const;
384
385 void gliftDebugShowALQ(DeferredLogger& deferred_logger);
386
388 GLiftProdWells& prod_wells,
389 GLiftOptWells& glift_wells,
390 GasLiftGroupInfo& group_info,
392 const int episodeIndex);
393
394 virtual void computePotentials(const std::size_t widx,
395 const WellState<Scalar>& well_state_copy,
396 std::string& exc_msg,
397 ExceptionType::ExcEnum& exc_type,
398 DeferredLogger& deferred_logger) = 0;
399
400 // Calculating well potentials for each well
401 void updateWellPotentials(const int reportStepIdx,
402 const bool onlyAfterEvent,
403 const SummaryConfig& summaryConfig,
404 DeferredLogger& deferred_logger);
405
407
408 void updateInjMult(DeferredLogger& deferred_logger);
409 void updateInjFCMult(DeferredLogger& deferred_logger);
410
411 void updateFiltrationParticleVolume(const double dt,
412 const std::size_t water_index);
413
414 // create the well container
415 virtual void createWellContainer(const int time_step) = 0;
416 virtual void initWellContainer(const int reportStepIdx) = 0;
417
418 virtual void calculateProductivityIndexValuesShutWells(const int reportStepIdx,
419 DeferredLogger& deferred_logger) = 0;
420 virtual void calculateProductivityIndexValues(DeferredLogger& deferred_logger) = 0;
421
422 void runWellPIScaling(const int reportStepIdx,
423 DeferredLogger& local_deferredLogger);
424
426 virtual int compressedIndexForInterior(int cartesian_cell_idx) const = 0;
427
428 std::vector<int> getCellsForConnections(const Well& well) const;
429 std::vector<std::vector<int>> getMaxWellConnections() const;
430
431 std::vector<std::string> getWellsForTesting(const int timeStepIdx,
432 const double simulationTime);
433
434 using WellTracerRates = std::map<std::pair<std::string, std::string>, Scalar>;
435 void assignWellTracerRates(data::Wells& wsrpt,
436 const WellTracerRates& wellTracerRates) const;
437
438 Schedule& schedule_;
439 const SummaryState& summaryState_;
440 const EclipseState& eclState_;
442
444 bool terminal_output_{false};
445 bool wells_active_{false};
446 bool network_active_{false};
449
450 std::optional<int> last_run_wellpi_{};
451
452 std::vector<Well> wells_ecl_;
453 std::vector<std::vector<PerforationData>> well_perf_data_;
454
457 {
458 public:
463 explicit ConnectionIndexMap(const std::size_t numConns)
464 : local_(numConns, -1)
465 {
466 this->global_.reserve(numConns);
467 this->open_.reserve(numConns);
468 }
469
477 void addActiveConnection(const int connIdx,
478 const bool connIsOpen)
479 {
480 this->local_[connIdx] =
481 static_cast<int>(this->global_.size());
482
483 this->global_.push_back(connIdx);
484
485 const auto open_conn_idx = connIsOpen
486 ? this->num_open_conns_++
487 : -1;
488
489 this->open_.push_back(open_conn_idx);
490 }
491
497 const std::vector<int>& local() const
498 {
499 return this->local_;
500 }
501
507 int global(const int connIdx) const
508 {
509 return this->global_[connIdx];
510 }
511
519 int open(const int connIdx) const
520 {
521 return this->open_[connIdx];
522 }
523
524 private:
528 std::vector<int> local_{};
529
532 std::vector<int> global_{};
533
535 std::vector<int> open_{};
536
538 int num_open_conns_{0};
539 };
540
541 std::vector<ConnectionIndexMap> conn_idx_map_{};
542 std::function<bool(const Well&)> not_on_process_{};
543
544 // a vector of all the wells.
545 std::vector<WellInterfaceGeneric*> well_container_generic_{};
546
547 std::vector<int> local_shut_wells_{};
548
549 std::vector<ParallelWellInfo> parallel_well_info_;
550 std::vector<std::reference_wrapper<ParallelWellInfo>> local_parallel_well_info_;
551
552 std::vector<WellProdIndexCalculator> prod_index_calc_;
554
555 std::vector<int> pvt_region_idx_;
556
557 mutable std::unordered_set<std::string> closed_this_step_;
558
559 GuideRate guideRate_;
560 std::unique_ptr<VFPProperties> vfp_properties_{};
561 std::map<std::string, Scalar> node_pressures_; // Storing network pressures for output.
562
563 // previous injection multiplier, it is used in the injection multiplier calculation for WINJMULT keyword
564 std::unordered_map<std::string, std::vector<Scalar>> prev_inj_multipliers_;
565
566 // Handling for filter cake injection multipliers
567 std::unordered_map<std::string, WellFilterCake<Scalar>> filter_cake_;
568
569 /*
570 The various wellState members should be accessed and modified
571 through the accessor functions wellState(), prevWellState(),
572 commitWellState(), resetWellState(), nupcolWellState() and
573 updateNupcolWellState().
574 */
578
579 bool glift_debug = false;
580
581 double last_glift_opt_time_ = -1.0;
582
584
585 // Store maps of group name and new group controls for output
586 std::map<std::string, std::string> switched_prod_groups_;
587 std::map<std::pair<std::string, Phase>, std::string> switched_inj_groups_;
588 // Store map of group name and close offending well for output
589 std::map<std::string, std::pair<std::string, std::string>> closed_offending_wells_;
590
591private:
592 WellInterfaceGeneric* getGenWell(const std::string& well_name);
593};
594
595
596} // namespace Opm
597
598#endif
Connection index mappings.
Definition: BlackoilWellModelGeneric.hpp:457
int open(const int connIdx) const
Definition: BlackoilWellModelGeneric.hpp:519
const std::vector< int > & local() const
Definition: BlackoilWellModelGeneric.hpp:497
int global(const int connIdx) const
Definition: BlackoilWellModelGeneric.hpp:507
ConnectionIndexMap(const std::size_t numConns)
Definition: BlackoilWellModelGeneric.hpp:463
void addActiveConnection(const int connIdx, const bool connIsOpen)
Definition: BlackoilWellModelGeneric.hpp:477
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:83
void gliftDebugShowALQ(DeferredLogger &deferred_logger)
std::unordered_map< std::string, WellFilterCake< Scalar > > filter_cake_
Definition: BlackoilWellModelGeneric.hpp:567
void updateNetworkActiveState(const int report_step)
Checks if network is active (at least one network well on prediction).
void updateNupcolWGState()
Definition: BlackoilWellModelGeneric.hpp:319
bool initial_step_
Definition: BlackoilWellModelGeneric.hpp:447
const Well & getWellEcl(const std::string &well_name) const
void assignNodeValues(std::map< std::string, data::NodeData > &nodevalues, const int reportStepIdx) const
void assignWellTracerRates(data::Wells &wsrpt, const WellTracerRates &wellTracerRates) const
const EclipseState & eclipseState() const
Definition: BlackoilWellModelGeneric.hpp:201
GuideRate guideRate_
Definition: BlackoilWellModelGeneric.hpp:559
std::map< std::pair< std::string, Phase >, std::string > switched_inj_groups_
Definition: BlackoilWellModelGeneric.hpp:587
Scalar updateNetworkPressures(const int reportStepIdx)
const SummaryState & summaryState_
Definition: BlackoilWellModelGeneric.hpp:439
const GuideRate & guideRate() const
Definition: BlackoilWellModelGeneric.hpp:205
bool reportStepStarts() const
Definition: BlackoilWellModelGeneric.hpp:207
virtual void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger)=0
void checkGconsaleLimits(const Group &group, WellState< Scalar > &well_state, const int reportStepIdx, DeferredLogger &deferred_logger)
std::vector< std::reference_wrapper< ParallelWellInfo > > createLocalParallelWellInfo(const std::vector< Well > &wells)
Create the parallel well information.
bool wellsActive() const
return true if wells are available in the reservoir
bool report_step_starts_
Definition: BlackoilWellModelGeneric.hpp:448
bool checkGroupHigherConstraints(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx)
GroupState< Scalar > & groupState()
Definition: BlackoilWellModelGeneric.hpp:145
std::vector< std::vector< int > > getMaxWellConnections() const
std::vector< WellInterfaceGeneric * > well_container_generic_
Definition: BlackoilWellModelGeneric.hpp:545
bool wellStructureChangedDynamically_
Definition: BlackoilWellModelGeneric.hpp:583
void updateAndCommunicateGroupData(const int reportStepIdx, const int iterationIdx)
bool networkActive() const
return true if network is active (at least one network well in prediction mode)
const Schedule & schedule() const
Definition: BlackoilWellModelGeneric.hpp:115
Schedule & schedule_
Definition: BlackoilWellModelGeneric.hpp:438
std::vector< Well > wells_ecl_
Definition: BlackoilWellModelGeneric.hpp:452
void updateWellPotentials(const int reportStepIdx, const bool onlyAfterEvent, const SummaryConfig &summaryConfig, DeferredLogger &deferred_logger)
bool forceShutWellByName(const std::string &wellname, const double simulation_time)
ParallelWBPCalculation wbpCalculationService_
Definition: BlackoilWellModelGeneric.hpp:553
std::map< std::string, Scalar > node_pressures_
Definition: BlackoilWellModelGeneric.hpp:561
std::vector< WellProdIndexCalculator > prod_index_calc_
Definition: BlackoilWellModelGeneric.hpp:552
void assignShutConnections(data::Wells &wsrpt, const int reportStepIndex) const
const WellTestState & wellTestState() const
Definition: BlackoilWellModelGeneric.hpp:149
void updateClosedWellsThisStep(const std::string &well_name) const
Definition: BlackoilWellModelGeneric.hpp:212
void assignGroupControl(const Group &group, data::GroupData &gdata) const
std::map< std::pair< std::string, std::string >, Scalar > WellTracerRates
Definition: BlackoilWellModelGeneric.hpp:434
void updateFiltrationParticleVolume(const double dt, const std::size_t water_index)
double last_glift_opt_time_
Definition: BlackoilWellModelGeneric.hpp:581
std::vector< int > local_shut_wells_
Definition: BlackoilWellModelGeneric.hpp:547
std::vector< int > pvt_region_idx_
Definition: BlackoilWellModelGeneric.hpp:555
void updateInjFCMult(DeferredLogger &deferred_logger)
virtual void createWellContainer(const int time_step)=0
void setWsolvent(const Group &group, const int reportStepIdx, Scalar wsolvent)
void updateEclWells(const int timeStepIdx, const SimulatorUpdate &sim_update, const SummaryState &st)
const WGState< Scalar > & prevWGState() const
Definition: BlackoilWellModelGeneric.hpp:289
std::vector< Well > getLocalWells(const int timeStepIdx) const
bool glift_debug
Definition: BlackoilWellModelGeneric.hpp:579
std::unordered_set< std::string > closed_this_step_
Definition: BlackoilWellModelGeneric.hpp:557
void gliftDebug(const std::string &msg, DeferredLogger &deferred_logger) const
data::GroupAndNetworkValues groupAndNetworkData(const int reportStepIdx) const
std::vector< std::reference_wrapper< ParallelWellInfo > > local_parallel_well_info_
Definition: BlackoilWellModelGeneric.hpp:550
virtual void initWellContainer(const int reportStepIdx)=0
bool wells_active_
Definition: BlackoilWellModelGeneric.hpp:445
bool shouldBalanceNetwork(const int reportStepIndex, const int iterationIdx) const
void initFromRestartFile(const RestartValue &restartValues, WellTestState wtestState, const std::size_t numCells, bool handle_ms_well)
const WellState< Scalar > & wellState() const
Definition: BlackoilWellModelGeneric.hpp:124
std::map< std::string, std::pair< std::string, std::string > > closed_offending_wells_
Definition: BlackoilWellModelGeneric.hpp:589
const std::vector< PerforationData > & perfData(const int well_idx) const
Definition: BlackoilWellModelGeneric.hpp:196
virtual void calcRates(const int fipnum, const int pvtreg, const std::vector< Scalar > &production_rates, std::vector< Scalar > &resv_coeff)=0
std::function< bool(const Well &)> not_on_process_
Definition: BlackoilWellModelGeneric.hpp:542
void prepareDeserialize(int report_step, const std::size_t numCells, bool handle_ms_well)
void commitWGState(WGState< Scalar > wgstate)
Definition: BlackoilWellModelGeneric.hpp:299
void commitWGState()
Definition: BlackoilWellModelGeneric.hpp:173
bool needPreStepNetworkRebalance(const int report_step) const
void assignGroupValues(const int reportStepIdx, std::map< std::string, data::GroupData > &gvalues) const
std::vector< std::vector< PerforationData > > well_perf_data_
Definition: BlackoilWellModelGeneric.hpp:453
std::map< std::string, const WellInterfaceGeneric * > GLiftProdWells
Definition: BlackoilWellModelGeneric.hpp:87
std::map< std::string, std::string > switched_prod_groups_
Definition: BlackoilWellModelGeneric.hpp:586
std::optional< int > last_run_wellpi_
Definition: BlackoilWellModelGeneric.hpp:450
bool hasTHPConstraints() const
Return true if any well has a THP constraint.
bool network_active_
Definition: BlackoilWellModelGeneric.hpp:446
std::unique_ptr< VFPProperties > vfp_properties_
Definition: BlackoilWellModelGeneric.hpp:560
const Parallel::Communication & comm() const
Definition: BlackoilWellModelGeneric.hpp:199
const WellState< Scalar > & prevWellState() const
Definition: BlackoilWellModelGeneric.hpp:284
const Parallel::Communication & comm_
Definition: BlackoilWellModelGeneric.hpp:441
WGState< Scalar > last_valid_wgstate_
Definition: BlackoilWellModelGeneric.hpp:576
bool wasDynamicallyShutThisTimeStep(const std::string &well_name) const
bool terminal_output_
Definition: BlackoilWellModelGeneric.hpp:444
std::vector< ParallelWellInfo > parallel_well_info_
Definition: BlackoilWellModelGeneric.hpp:549
WGState< Scalar > nupcol_wgstate_
Definition: BlackoilWellModelGeneric.hpp:577
std::vector< const WellInterfaceGeneric * > genericWells() const
Definition: BlackoilWellModelGeneric.hpp:118
std::vector< ConnectionIndexMap > conn_idx_map_
Definition: BlackoilWellModelGeneric.hpp:541
const GroupState< Scalar > & groupState() const
Definition: BlackoilWellModelGeneric.hpp:117
virtual void calculateProductivityIndexValues(DeferredLogger &deferred_logger)=0
bool operator==(const BlackoilWellModelGeneric &rhs) const
Definition: BlackoilWellModelGeneric.hpp:238
WGState< Scalar > active_wgstate_
Definition: BlackoilWellModelGeneric.hpp:575
BlackoilWellModelGeneric(Schedule &schedule, const SummaryState &summaryState, const EclipseState &eclState, const PhaseUsage &phase_usage, const Parallel::Communication &comm)
WellState< Scalar > & wellState()
Definition: BlackoilWellModelGeneric.hpp:132
WellTestState & wellTestState()
Definition: BlackoilWellModelGeneric.hpp:147
std::map< std::string, std::unique_ptr< GasLiftSingleWellGeneric > > GLiftOptWells
Definition: BlackoilWellModelGeneric.hpp:86
std::vector< int > getCellsForConnections(const Well &well) const
const PhaseUsage & phaseUsage() const
Definition: BlackoilWellModelGeneric.hpp:116
void updateWsolvent(const Group &group, const int reportStepIdx, const WellState< Scalar > &wellState)
virtual ~BlackoilWellModelGeneric()=default
void calculateEfficiencyFactors(const int reportStepIdx)
virtual void computePotentials(const std::size_t widx, const WellState< Scalar > &well_state_copy, std::string &exc_msg, ExceptionType::ExcEnum &exc_type, DeferredLogger &deferred_logger)=0
const SummaryState & summaryState() const
Definition: BlackoilWellModelGeneric.hpp:203
Scalar wellPI(const std::string &well_name) const
void gasLiftOptimizationStage2(DeferredLogger &deferred_logger, GLiftProdWells &prod_wells, GLiftOptWells &glift_wells, GasLiftGroupInfo &group_info, GLiftWellStateMap &map, const int episodeIndex)
std::unordered_map< std::string, std::vector< Scalar > > prev_inj_multipliers_
Definition: BlackoilWellModelGeneric.hpp:564
void serializeOp(Serializer &serializer)
Definition: BlackoilWellModelGeneric.hpp:219
const EclipseState & eclState_
Definition: BlackoilWellModelGeneric.hpp:440
std::map< std::string, std::unique_ptr< GasLiftWellState > > GLiftWellStateMap
Definition: BlackoilWellModelGeneric.hpp:88
bool wasDynamicallyShutThisTimeStep(const int well_index) const
void checkGEconLimits(const Group &group, const double simulation_time, const int report_step_idx, DeferredLogger &deferred_logger)
void runWellPIScaling(const int reportStepIdx, DeferredLogger &local_deferredLogger)
bool hasWell(const std::string &wname) const
void resetWGState()
Definition: BlackoilWellModelGeneric.hpp:309
Scalar wellPI(const int well_index) const
PhaseUsage phase_usage_
Definition: BlackoilWellModelGeneric.hpp:443
const WellState< Scalar > & nupcolWellState() const
Definition: BlackoilWellModelGeneric.hpp:141
virtual void calcInjRates(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff)=0
void updateInjMult(DeferredLogger &deferred_logger)
std::vector< std::string > getWellsForTesting(const int timeStepIdx, const double simulationTime)
virtual int compressedIndexForInterior(int cartesian_cell_idx) const =0
get compressed index for interior cells (-1, otherwise
Definition: DeferredLogger.hpp:57
Definition: GasLiftGroupInfo.hpp:45
Definition: GroupState.hpp:35
Definition: ParallelWBPCalculation.hpp:50
Definition: WellInterfaceGeneric.hpp:50
Definition: WellState.hpp:62
ExcEnum
Definition: DeferredLogger.hpp:45
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: BlackoilPhases.hpp:27
Definition: BlackoilPhases.hpp:46
Definition: WGState.hpp:39