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
46
47#include <algorithm>
48#include <cstddef>
49#include <functional>
50#include <map>
51#include <memory>
52#include <optional>
53#include <string>
54#include <unordered_map>
55#include <unordered_set>
56#include <vector>
57
58
59namespace Opm {
60 class DeferredLogger;
61 class EclipseState;
62 template<typename Scalar, typename IndexTraits> class BlackoilWellModelGasLiftGeneric;
63 template<typename Scalar, typename IndexTraits> class BlackoilWellModelNetworkGeneric;
64 template<typename Scalar, typename IndexTraits> class GasLiftGroupInfo;
65 template<typename Scalar, typename IndexTraits> class GasLiftSingleWellGeneric;
66 template<class Scalar> class GasLiftWellState;
67 class Group;
68 class GuideRateConfig;
69 class RestartValue;
70 class Schedule;
71 struct SimulatorUpdate;
72 class SummaryConfig;
73 template<typename Scalar, typename IndexTraits> class VFPProperties;
74 template<typename Scalar, typename IndexTraits> class WellInterfaceGeneric;
75 template<typename Scalar, typename IndexTraits> class WellState;
76} // namespace Opm
77
78namespace Opm { namespace data {
79 struct GroupData;
80 struct GroupGuideRates;
81 class GroupAndNetworkValues;
82 struct NodeData;
83}} // namespace Opm::data
84
85namespace Opm::Parameters {
86
87struct EnableTerminalOutput { static constexpr bool value = true; };
88
89} // namespace Opm::Parameters
90
91namespace Opm {
92
94template<typename Scalar, typename IndexTraits>
96{
98public:
102 const SummaryState& summaryState,
103 const EclipseState& eclState,
104 const PhaseUsageInfo<IndexTraits>& phase_usage,
106
107 virtual ~BlackoilWellModelGeneric() = default;
108 virtual int compressedIndexForInteriorLGR([[maybe_unused]] const std::string& lgr_tag,
109 [[maybe_unused]] const Connection& conn) const
110 {
111 throw std::runtime_error("compressedIndexForInteriorLGR not implemented");
112 }
113
114 int numLocalWells() const;
115 int numLocalWellsEnd() const;
117 int numPhases() const;
118
120 bool wellsActive() const;
121
123 bool hasLocalWell(const std::string& wname) const;
124
126 bool hasOpenLocalWell(const std::string& well_name) const;
127
128 // whether there exists any multisegment well open on this process
129 bool anyMSWellOpenLocal() const;
130
131 const std::vector<Well>& eclWells() const
132 { return wells_ecl_; }
133
134 bool terminalOutput() const
135 { return terminal_output_; }
136
137 const Well& getWellEcl(const std::string& well_name) const;
138 std::vector<Well> getLocalWells(const int timeStepIdx) const;
139 const Schedule& schedule() const { return schedule_; }
141 const GroupState<Scalar>& groupState() const { return this->active_wgstate_.group_state; }
142 std::vector<const WellInterfaceGeneric<Scalar, IndexTraits>*> genericWells() const
143 { return {well_container_generic_.begin(), well_container_generic_.end()}; }
144
145 std::vector<WellInterfaceGeneric<Scalar, IndexTraits>*> genericWells()
146 { return well_container_generic_; }
147
148 /*
149 Immutable version of the currently active wellstate.
150 */
152 {
153 return this->active_wgstate_.well_state;
154 }
155
156 /*
157 Mutable version of the currently active wellstate.
158 */
160 {
161 return this->active_wgstate_.well_state;
162 }
163
164 /*
165 Will return the currently active nupcolWellState; must update
166 the internal nupcol wellstate with updateNupcolWGState() first.
167
168 Both const and non-const accessors are provided. The non-const
169 accessor is required for the WellStateGuard pattern and pushWellState()
170 in WellGroupHelper, which temporarily switches WellGroupHelper to use this state.
171 */
173 {
174 return this->nupcol_wgstate_.well_state;
175 }
177 {
178 return this->nupcol_wgstate_.well_state;
179 }
180 GroupState<Scalar>& groupState() { return this->active_wgstate_.group_state; }
181
182 WellTestState& wellTestState() { return this->active_wgstate_.well_test_state; }
183
184 const WellTestState& wellTestState() const { return this->active_wgstate_.well_test_state; }
185
186 Scalar wellPI(const int well_index) const;
187 Scalar wellPI(const std::string& well_name) const;
188
189 void updateEclWells(const int timeStepIdx,
190 const SimulatorUpdate& sim_update,
191 const SummaryState& st);
192
193 void initFromRestartFile(const RestartValue& restartValues,
194 std::unique_ptr<WellTestState> wtestState,
195 const std::size_t numCells,
196 bool handle_ms_well,
197 bool enable_distributed_wells);
198
199 void prepareDeserialize(int report_step,
200 const std::size_t numCells,
201 bool handle_ms_well,
202 bool enable_distributed_wells);
203
204 /*
205 Will assign the internal member last_valid_well_state_ to the
206 current value of the this->active_well_state_. The state stored
207 with storeWellState() can then subsequently be recovered with the
208 resetWellState() method.
209 */
211
212 data::GroupAndNetworkValues groupAndNetworkData(const int reportStepIdx) const;
213
216 bool forceShutWellByName(const std::string& wellname,
217 const double simulation_time,
218 const bool dont_shut_grup_wells);
219
220 const std::vector<PerforationData<Scalar>>& perfData(const int well_idx) const
221 { return well_perf_data_[well_idx]; }
222
223 const Parallel::Communication& comm() const { return comm_; }
224
225 const EclipseState& eclipseState() const { return eclState_; }
226
227 const SummaryState& summaryState() const { return summaryState_; }
228
229 const GuideRate& guideRate() const { return guideRate_; }
230 GuideRate& guideRate() { return guideRate_; }
231
232 const std::map<std::string, double>& wellOpenTimes() const { return well_open_times_; }
233 const std::map<std::string, double>& wellCloseTimes() const { return well_close_times_; }
234 const WellGroupEvents& reportStepStartEvents() const { return report_step_start_events_; }
235
236 std::vector<int> getCellsForConnections(const Well& well) const;
237
238 bool reportStepStarts() const { return report_step_starts_; }
239
240 void updateClosedWellsThisStep(const std::string& well_name) const
241 {
242 this->closed_this_step_.insert(well_name);
243 }
244 bool wasDynamicallyShutThisTimeStep(const std::string& well_name) const;
245
246 void logPrimaryVars() const;
247
248 template<class Serializer>
249 void serializeOp(Serializer& serializer)
250 {
251 serializer(initial_step_);
252 serializer(report_step_starts_);
253 serializer(last_run_wellpi_);
254 serializer(local_shut_wells_);
255 serializer(closed_this_step_);
256 serializer(guideRate_);
257 serializer(genNetwork_);
258 serializer(prev_inj_multipliers_);
259 serializer(active_wgstate_);
260 serializer(last_valid_wgstate_);
261 serializer(nupcol_wgstate_);
262 serializer(switched_prod_groups_);
263 serializer(switched_inj_groups_);
264 serializer(closed_offending_wells_);
265 serializer(gen_gaslift_);
266 }
267
268 bool operator==(const BlackoilWellModelGeneric& rhs) const;
269
271 parallelWellInfo(const std::size_t idx) const
272 { return local_parallel_well_info_[idx].get(); }
273
274 bool isOwner(const std::string& wname) const
275 {
276 return this->parallelWellSatisfies
277 (wname, [](const auto& pwInfo) { return pwInfo.isOwner(); });
278 }
279
280 bool hasLocalCells(const std::string& wname) const
281 {
282 return this->parallelWellSatisfies
283 (wname, [](const auto& pwInfo) { return pwInfo.hasLocalCells(); });
284 }
285
286 const ConnectionIndexMap& connectionIndexMap(const std::size_t idx)
287 { return conn_idx_map_[idx]; }
288
291
293 {
294 return *vfp_properties_;
295 }
296
297 void updateAndCommunicateGroupData(const int reportStepIdx,
298 const int iterationIdx,
299 const Scalar tol_nupcol,
300 // we only want to update the wellgroup target
301 // after the groups have found their controls
302 const bool update_wellgrouptarget,
303 DeferredLogger& deferred_logger);
304
305 const EclipseState& eclState() const
306 { return eclState_; }
307
308protected:
309 /*
310 The dynamic state of the well model is maintained with an instance
311 of the WellState class. Currently we have
312 three different wellstate instances:
313
314 1. The currently active wellstate is in the active_well_state_
315 member. That is the state which is mutated by the simulator.
316
317 2. In the case timestep fails to converge and we must go back and
318 try again with a smaller timestep we need to recover the last
319 valid wellstate. This is maintained with the
320 last_valid_well_state_ member and the functions
321 commitWGState() and resetWellState().
322
323 3. For the NUPCOL functionality we should either use the
324 currently active wellstate or a wellstate frozen at max
325 nupcol iterations. This is handled with the member
326 nupcol_well_state_ and the updateNupcolWGState() function.
327 */
328
329 /*
330 Will return the last good wellstate. This is typcially used when
331 initializing a new report step where the Schedule object might
332 have introduced new wells. The wellstate returned by
333 prevWellState() must have been stored with the commitWGState()
334 function first.
335 */
337 {
338 return this->last_valid_wgstate_.well_state;
339 }
340
342 {
343 return this->last_valid_wgstate_;
344 }
345
346 /*
347 Will store a copy of the input argument well_state in the
348 last_valid_well_state_ member, that state can then be recovered
349 with a subsequent call to resetWellState().
350 */
352 {
353 this->last_valid_wgstate_ = std::move(wgstate);
354 }
355
356 /*
357 Will update the internal variable active_well_state_ to whatever
358 was stored in the last_valid_well_state_ member. This function
359 works in pair with commitWellState() which should be called first.
360 */
362 {
364 this->genNetwork_.resetState();
365 // Update helper pointers to reference the restored active state
366 this->group_state_helper_.updateState(this->wellState(), this->groupState());
367 }
368
369 /*
370 Will store the current active wellstate in the nupcol_well_state_
371 member. This can then be subsequently retrieved with accessor
372 nupcolWellState().
373 */
375 {
376 this->nupcol_wgstate_ = this->active_wgstate_;
377 }
378
379 void reportGroupSwitching(DeferredLogger& local_deferredLogger) const;
380
383 std::vector<std::reference_wrapper<ParallelWellInfo<Scalar>>>
384 createLocalParallelWellInfo(const std::vector<Well>& wells);
385
388
389 bool wasDynamicallyShutThisTimeStep(const int well_index) const;
390
391 void updateWsolvent(const Group& group,
392 const int reportStepIdx,
394 void setWsolvent(const Group& group,
395 const int reportStepIdx,
396 Scalar wsolvent);
397 virtual void calcResvCoeff(const int fipnum,
398 const int pvtreg,
399 const std::vector<Scalar>& production_rates,
400 std::vector<Scalar>& resv_coeff) = 0;
401 virtual void calcInjResvCoeff(const int fipnum,
402 const int pvtreg,
403 std::vector<Scalar>& resv_coeff) = 0;
404
409 void assignDynamicWellStatus(data::Wells& wsrpt) const;
410
422 void assignShutConnections(data::Wells& wsrpt,
423 const int reportStepIndex) const;
424
425 void assignWellTargets(data::Wells& wsrpt) const;
426 void assignProductionWellTargets(const Well& well, data::WellControlLimits& limits) const;
427 void assignInjectionWellTargets(const Well& well, data::WellControlLimits& limits) const;
428
429 void assignGroupControl(const Group& group,
430 data::GroupData& gdata) const;
431 void assignGroupValues(const int reportStepIdx,
432 std::map<std::string, data::GroupData>& gvalues) const;
433
434 void calculateEfficiencyFactors(const int reportStepIdx);
435
436 void checkGconsaleLimits(const Group& group,
438 const int reportStepIdx,
439 DeferredLogger& deferred_logger);
440
441 void checkGEconLimits(const Group& group,
442 const double simulation_time,
443 const int report_step_idx,
444 DeferredLogger& deferred_logger);
445
446 bool checkGroupHigherConstraints(const Group& group,
447 DeferredLogger& deferred_logger,
448 const int reportStepIdx,
449 const int max_number_of_group_switch,
450 const bool update_group_switching_log);
451
453
455
456 virtual void computePotentials(const std::size_t widx,
457 const WellState<Scalar, IndexTraits>& well_state_copy,
458 std::string& exc_msg,
459 ExceptionType::ExcEnum& exc_type,
460 DeferredLogger& deferred_logger) = 0;
461
462 // Calculating well potentials for each well
463 void updateWellPotentials(const int reportStepIdx,
464 const bool onlyAfterEvent,
465 const SummaryConfig& summaryConfig,
466 DeferredLogger& deferred_logger);
467
469
470 void updateInjMult(DeferredLogger& deferred_logger);
471 void updateInjFCMult(DeferredLogger& deferred_logger);
472
473 void updateFiltrationModelsPostStep(const double dt,
474 const std::size_t water_index,
475 DeferredLogger& deferred_logger);
476
478
479 // create the well container
480 virtual void createWellContainer(const int time_step) = 0;
481 virtual void initWellContainer(const int reportStepIdx) = 0;
482
483 virtual void calculateProductivityIndexValuesShutWells(const int reportStepIdx,
484 DeferredLogger& deferred_logger) = 0;
485 virtual void calculateProductivityIndexValues(DeferredLogger& deferred_logger) = 0;
486
487 void runWellPIScaling(const int reportStepIdx,
488 DeferredLogger& local_deferredLogger);
489
491 virtual int compressedIndexForInterior(int cartesian_cell_idx) const = 0;
492
493 std::vector<std::vector<int>> getMaxWellConnections() const;
494
495 std::vector<std::string> getWellsForTesting(const int timeStepIdx,
496 const double simulationTime);
497
498 using WellTracerRates = std::unordered_map<int, std::vector<WellTracerRate<Scalar>>>;
499 void assignWellTracerRates(data::Wells& wsrpt,
500 const WellTracerRates& wellTracerRates,
501 const unsigned reportStep) const;
502
503 using MswTracerRates = std::unordered_map<int, std::vector<MSWellTracerRate<Scalar>>>;
504 void assignMswTracerRates(data::Wells& wsrpt,
505 const MswTracerRates& mswTracerRates,
506 const unsigned reportStep) const;
507
508 void assignMassGasRate(data::Wells& wsrpt,
509 const Scalar gasDensity) const;
510
511 Schedule& schedule_;
512
513 const SummaryState& summaryState_;
514 const EclipseState& eclState_;
518
520 bool terminal_output_{false};
521 bool wells_active_{false};
524
525 std::optional<int> last_run_wellpi_{};
526
527 std::vector<Well> wells_ecl_;
528 std::vector<std::vector<PerforationData<Scalar>>> well_perf_data_;
529
530 // Times at which wells were opened (for WCYCLE)
531 std::map<std::string, double> well_open_times_;
532
533 // Times at which wells were shut (for WCYCLE)
534 std::map<std::string, double> well_close_times_;
535
536 std::vector<ConnectionIndexMap> conn_idx_map_{};
537 std::function<bool(const std::string&)> not_on_process_{};
538
539 // a vector of all the wells.
540 std::vector<WellInterfaceGeneric<Scalar, IndexTraits>*> well_container_generic_{};
541
542 std::vector<int> local_shut_wells_{};
543
544 std::vector<ParallelWellInfo<Scalar>> parallel_well_info_;
545 std::vector<std::reference_wrapper<ParallelWellInfo<Scalar>>> local_parallel_well_info_;
546
547 std::vector<WellProdIndexCalculator<Scalar>> prod_index_calc_;
548
549 std::vector<int> pvt_region_idx_;
550
551 mutable std::unordered_set<std::string> closed_this_step_;
552
553 GuideRate guideRate_;
554 std::unique_ptr<VFPProperties<Scalar, IndexTraits>> vfp_properties_{};
555
556 // previous injection multiplier, it is used in the injection multiplier calculation for WINJMULT keyword
557 std::unordered_map<std::string, std::vector<Scalar>> prev_inj_multipliers_;
558
559 // Handling for filter cake injection multipliers
560 std::unordered_map<std::string, WellFilterCake<Scalar, IndexTraits>> filter_cake_;
561
562 /*
563 The various wellState members should be accessed and modified
564 through the accessor functions wellState(), prevWellState(),
565 commitWellState(), resetWellState(), nupcolWellState() and
566 updateNupcolWGState().
567 */
573
575
576 // Store maps of group name and new group controls for output
577 std::map<std::string, std::vector<Group::ProductionCMode>> switched_prod_groups_;
578 std::map<std::string, std::array<std::vector<Group::InjectionCMode>, 3>> switched_inj_groups_;
579 // Store map of group name and close offending well for output
580 std::map<std::string, std::pair<std::string, std::string>> closed_offending_wells_;
581
583
584private:
585 WellInterfaceGeneric<Scalar, IndexTraits>* getGenWell(const std::string& well_name);
586
587 template <typename Iter, typename Body>
588 void wellUpdateLoop(Iter first, Iter last, const int timeStepIdx, Body&& body);
589
590 void updateEclWellsConstraints(const int timeStepIdx,
591 const SimulatorUpdate& sim_update,
592 const SummaryState& st);
593
594 void updateEclWellsCTFFromAction(const int timeStepIdx,
595 const SimulatorUpdate& sim_update);
596
610 template <typename Predicate>
611 bool parallelWellSatisfies(const std::string& wname, Predicate&& p) const
612 {
613 auto pwInfoPos = std::find_if(this->parallel_well_info_.begin(),
614 this->parallel_well_info_.end(),
615 [&wname](const auto& pwInfo)
616 { return pwInfo.name() == wname; });
617
618 return (pwInfoPos != this->parallel_well_info_.end())
619 && p(*pwInfoPos);
620 }
621
637 template <typename LoopBody>
638 void loopOwnedWells(LoopBody&& loopBody) const;
639};
640
641} // namespace Opm
642
643#endif
Definition: BlackoilWellModelGasLift.hpp:42
Class for handling the blackoil well model.
Definition: BlackoilWellModelGeneric.hpp:96
std::map< std::string, std::vector< Group::ProductionCMode > > switched_prod_groups_
Definition: BlackoilWellModelGeneric.hpp:577
void updateAndCommunicateGroupData(const int reportStepIdx, const int iterationIdx, const Scalar tol_nupcol, const bool update_wellgrouptarget, DeferredLogger &deferred_logger)
GroupStateHelperType & groupStateHelper()
Definition: BlackoilWellModelGeneric.hpp:289
const EclipseState & eclState_
Definition: BlackoilWellModelGeneric.hpp:514
const PhaseUsageInfo< IndexTraits > & phaseUsage() const
Definition: BlackoilWellModelGeneric.hpp:140
const WGState< Scalar, IndexTraits > & prevWGState() const
Definition: BlackoilWellModelGeneric.hpp:341
std::vector< Well > getLocalWells(const int timeStepIdx) const
const ConnectionIndexMap & connectionIndexMap(const std::size_t idx)
Definition: BlackoilWellModelGeneric.hpp:286
WellState< Scalar, IndexTraits > & nupcolWellState()
Definition: BlackoilWellModelGeneric.hpp:176
const GroupState< Scalar > & groupState() const
Definition: BlackoilWellModelGeneric.hpp:141
const PhaseUsageInfo< IndexTraits > & phase_usage_info_
Definition: BlackoilWellModelGeneric.hpp:519
void assignDynamicWellStatus(data::Wells &wsrpt) const
GroupStateHelperType group_state_helper_
Definition: BlackoilWellModelGeneric.hpp:571
void serializeOp(Serializer &serializer)
Definition: BlackoilWellModelGeneric.hpp:249
bool reportStepStarts() const
Definition: BlackoilWellModelGeneric.hpp:238
bool operator==(const BlackoilWellModelGeneric &rhs) const
WGState< Scalar, IndexTraits > last_valid_wgstate_
Definition: BlackoilWellModelGeneric.hpp:569
const Parallel::Communication & comm() const
Definition: BlackoilWellModelGeneric.hpp:223
virtual int compressedIndexForInteriorLGR(const std::string &lgr_tag, const Connection &conn) const
Definition: BlackoilWellModelGeneric.hpp:108
GuideRate guideRate_
Definition: BlackoilWellModelGeneric.hpp:553
virtual int compressedIndexForInterior(int cartesian_cell_idx) const =0
get compressed index for interior cells (-1, otherwise
std::vector< WellInterfaceGeneric< Scalar, IndexTraits > * > well_container_generic_
Definition: BlackoilWellModelGeneric.hpp:540
WellGroupEvents report_step_start_events_
Well group events at start of report step.
Definition: BlackoilWellModelGeneric.hpp:572
void resetWGState()
Definition: BlackoilWellModelGeneric.hpp:361
void assignWellTargets(data::Wells &wsrpt) const
virtual void calcResvCoeff(const int fipnum, const int pvtreg, const std::vector< Scalar > &production_rates, std::vector< Scalar > &resv_coeff)=0
WGState< Scalar, IndexTraits > active_wgstate_
Definition: BlackoilWellModelGeneric.hpp:568
Scalar wellPI(const int well_index) const
virtual ~BlackoilWellModelGeneric()=default
std::vector< Well > wells_ecl_
Definition: BlackoilWellModelGeneric.hpp:527
bool hasLocalWell(const std::string &wname) const
Returns true if well is defined and has connections on current rank.
void assignInjectionWellTargets(const Well &well, data::WellControlLimits &limits) const
const ParallelWellInfo< Scalar > & parallelWellInfo(const std::size_t idx) const
Definition: BlackoilWellModelGeneric.hpp:271
BlackoilWellModelWBP< Scalar, IndexTraits > wbp_
Definition: BlackoilWellModelGeneric.hpp:517
const WellGroupEvents & reportStepStartEvents() const
Definition: BlackoilWellModelGeneric.hpp:234
std::vector< const WellInterfaceGeneric< Scalar, IndexTraits > * > genericWells() const
Definition: BlackoilWellModelGeneric.hpp:142
const std::map< std::string, double > & wellCloseTimes() const
Definition: BlackoilWellModelGeneric.hpp:233
const std::vector< Well > & eclWells() const
Definition: BlackoilWellModelGeneric.hpp:131
void updateNupcolWGState()
Definition: BlackoilWellModelGeneric.hpp:374
void updateFiltrationModelsPostStep(const double dt, const std::size_t water_index, DeferredLogger &deferred_logger)
const EclipseState & eclipseState() const
Definition: BlackoilWellModelGeneric.hpp:225
bool wells_active_
Definition: BlackoilWellModelGeneric.hpp:521
WGState< Scalar, IndexTraits > nupcol_wgstate_
Definition: BlackoilWellModelGeneric.hpp:570
std::unordered_map< int, std::vector< MSWellTracerRate< Scalar > > > MswTracerRates
Definition: BlackoilWellModelGeneric.hpp:503
std::unique_ptr< VFPProperties< Scalar, IndexTraits > > vfp_properties_
Definition: BlackoilWellModelGeneric.hpp:554
bool report_step_starts_
Definition: BlackoilWellModelGeneric.hpp:523
std::vector< std::reference_wrapper< ParallelWellInfo< Scalar > > > local_parallel_well_info_
Definition: BlackoilWellModelGeneric.hpp:545
const Well & getWellEcl(const std::string &well_name) const
const Parallel::Communication & comm_
Definition: BlackoilWellModelGeneric.hpp:515
std::vector< WellInterfaceGeneric< Scalar, IndexTraits > * > genericWells()
Definition: BlackoilWellModelGeneric.hpp:145
std::vector< int > getCellsForConnections(const Well &well) const
void updateWsolvent(const Group &group, const int reportStepIdx, const WellState< Scalar, IndexTraits > &wellState)
std::unordered_map< int, std::vector< WellTracerRate< Scalar > > > WellTracerRates
Definition: BlackoilWellModelGeneric.hpp:498
WellTestState & wellTestState()
Definition: BlackoilWellModelGeneric.hpp:182
bool wellStructureChangedDynamically_
Definition: BlackoilWellModelGeneric.hpp:574
std::vector< std::vector< PerforationData< Scalar > > > well_perf_data_
Definition: BlackoilWellModelGeneric.hpp:528
const GroupStateHelperType & groupStateHelper() const
Definition: BlackoilWellModelGeneric.hpp:290
bool wellsActive() const
return true if wells are available in the reservoir
bool forceShutWellByName(const std::string &wellname, const double simulation_time, const bool dont_shut_grup_wells)
bool terminal_output_
Definition: BlackoilWellModelGeneric.hpp:520
const WellTestState & wellTestState() const
Definition: BlackoilWellModelGeneric.hpp:184
void reportGroupSwitching(DeferredLogger &local_deferredLogger) const
std::unordered_map< std::string, WellFilterCake< Scalar, IndexTraits > > filter_cake_
Definition: BlackoilWellModelGeneric.hpp:560
std::map< std::string, double > well_close_times_
Definition: BlackoilWellModelGeneric.hpp:534
virtual void calculateProductivityIndexValues(DeferredLogger &deferred_logger)=0
void checkGconsaleLimits(const Group &group, WellState< Scalar, IndexTraits > &well_state, const int reportStepIdx, DeferredLogger &deferred_logger)
void commitWGState(WGState< Scalar, IndexTraits > wgstate)
Definition: BlackoilWellModelGeneric.hpp:351
data::GroupAndNetworkValues groupAndNetworkData(const int reportStepIdx) const
void assignShutConnections(data::Wells &wsrpt, const int reportStepIndex) const
std::vector< int > pvt_region_idx_
Definition: BlackoilWellModelGeneric.hpp:549
void setWsolvent(const Group &group, const int reportStepIdx, Scalar wsolvent)
virtual void createWellContainer(const int time_step)=0
std::function< bool(const std::string &)> not_on_process_
Definition: BlackoilWellModelGeneric.hpp:537
const SummaryState & summaryState() const
Definition: BlackoilWellModelGeneric.hpp:227
bool wasDynamicallyShutThisTimeStep(const std::string &well_name) const
virtual void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger &deferred_logger)=0
const GuideRate & guideRate() const
Definition: BlackoilWellModelGeneric.hpp:229
bool checkGroupHigherConstraints(const Group &group, DeferredLogger &deferred_logger, const int reportStepIdx, const int max_number_of_group_switch, const bool update_group_switching_log)
std::optional< int > last_run_wellpi_
Definition: BlackoilWellModelGeneric.hpp:525
const SummaryState & summaryState_
Definition: BlackoilWellModelGeneric.hpp:513
std::map< std::string, std::array< std::vector< Group::InjectionCMode >, 3 > > switched_inj_groups_
Definition: BlackoilWellModelGeneric.hpp:578
std::map< std::string, std::pair< std::string, std::string > > closed_offending_wells_
Definition: BlackoilWellModelGeneric.hpp:580
void updateClosedWellsThisStep(const std::string &well_name) const
Definition: BlackoilWellModelGeneric.hpp:240
const WellState< Scalar, IndexTraits > & wellState() const
Definition: BlackoilWellModelGeneric.hpp:151
std::vector< int > local_shut_wells_
Definition: BlackoilWellModelGeneric.hpp:542
void updateFiltrationModelsPreStep(DeferredLogger &deferred_logger)
void calculateEfficiencyFactors(const int reportStepIdx)
virtual void calcInjResvCoeff(const int fipnum, const int pvtreg, std::vector< Scalar > &resv_coeff)=0
bool hasLocalCells(const std::string &wname) const
Definition: BlackoilWellModelGeneric.hpp:280
void assignMswTracerRates(data::Wells &wsrpt, const MswTracerRates &mswTracerRates, const unsigned reportStep) const
void updateInjFCMult(DeferredLogger &deferred_logger)
std::vector< ParallelWellInfo< Scalar > > parallel_well_info_
Definition: BlackoilWellModelGeneric.hpp:544
void assignProductionWellTargets(const Well &well, data::WellControlLimits &limits) const
void checkGEconLimits(const Group &group, const double simulation_time, const int report_step_idx, DeferredLogger &deferred_logger)
void assignGroupValues(const int reportStepIdx, std::map< std::string, data::GroupData > &gvalues) const
const Schedule & schedule() const
Definition: BlackoilWellModelGeneric.hpp:139
WellState< Scalar, IndexTraits > & wellState()
Definition: BlackoilWellModelGeneric.hpp:159
BlackoilWellModelNetworkGeneric< Scalar, IndexTraits > & genNetwork_
Definition: BlackoilWellModelGeneric.hpp:582
BlackoilWellModelGeneric(Schedule &schedule, BlackoilWellModelGasLiftGeneric< Scalar, IndexTraits > &gaslift, BlackoilWellModelNetworkGeneric< Scalar, IndexTraits > &network, const SummaryState &summaryState, const EclipseState &eclState, const PhaseUsageInfo< IndexTraits > &phase_usage, const Parallel::Communication &comm)
void assignGroupControl(const Group &group, data::GroupData &gdata) const
void assignMassGasRate(data::Wells &wsrpt, const Scalar gasDensity) const
bool hasOpenLocalWell(const std::string &well_name) const
Returns true if well is defined, open and has connections on current rank.
GuideRate & guideRate()
Definition: BlackoilWellModelGeneric.hpp:230
const EclipseState & eclState() const
Definition: BlackoilWellModelGeneric.hpp:305
BlackoilWellModelGasLiftGeneric< Scalar, IndexTraits > & gen_gaslift_
Definition: BlackoilWellModelGeneric.hpp:516
const std::map< std::string, double > & wellOpenTimes() const
Definition: BlackoilWellModelGeneric.hpp:232
void initFromRestartFile(const RestartValue &restartValues, std::unique_ptr< WellTestState > wtestState, const std::size_t numCells, bool handle_ms_well, bool enable_distributed_wells)
const WellState< Scalar, IndexTraits > & prevWellState() const
Definition: BlackoilWellModelGeneric.hpp:336
virtual void computePotentials(const std::size_t widx, const WellState< Scalar, IndexTraits > &well_state_copy, std::string &exc_msg, ExceptionType::ExcEnum &exc_type, DeferredLogger &deferred_logger)=0
const VFPProperties< Scalar, IndexTraits > & getVFPProperties() const
Definition: BlackoilWellModelGeneric.hpp:292
std::vector< std::vector< int > > getMaxWellConnections() const
void assignWellTracerRates(data::Wells &wsrpt, const WellTracerRates &wellTracerRates, const unsigned reportStep) const
std::vector< std::string > getWellsForTesting(const int timeStepIdx, const double simulationTime)
void updateEclWells(const int timeStepIdx, const SimulatorUpdate &sim_update, const SummaryState &st)
void updateInjMult(DeferredLogger &deferred_logger)
bool initial_step_
Definition: BlackoilWellModelGeneric.hpp:522
bool terminalOutput() const
Definition: BlackoilWellModelGeneric.hpp:134
void prepareDeserialize(int report_step, const std::size_t numCells, bool handle_ms_well, bool enable_distributed_wells)
GroupState< Scalar > & groupState()
Definition: BlackoilWellModelGeneric.hpp:180
bool wasDynamicallyShutThisTimeStep(const int well_index) const
const std::vector< PerforationData< Scalar > > & perfData(const int well_idx) const
Definition: BlackoilWellModelGeneric.hpp:220
std::vector< std::reference_wrapper< ParallelWellInfo< Scalar > > > createLocalParallelWellInfo(const std::vector< Well > &wells)
Create the parallel well information.
std::unordered_map< std::string, std::vector< Scalar > > prev_inj_multipliers_
Definition: BlackoilWellModelGeneric.hpp:557
Schedule & schedule_
Definition: BlackoilWellModelGeneric.hpp:511
const WellState< Scalar, IndexTraits > & nupcolWellState() const
Definition: BlackoilWellModelGeneric.hpp:172
std::vector< ConnectionIndexMap > conn_idx_map_
Definition: BlackoilWellModelGeneric.hpp:536
virtual void initWellContainer(const int reportStepIdx)=0
Scalar wellPI(const std::string &well_name) const
void updateWellPotentials(const int reportStepIdx, const bool onlyAfterEvent, const SummaryConfig &summaryConfig, DeferredLogger &deferred_logger)
std::map< std::string, double > well_open_times_
Definition: BlackoilWellModelGeneric.hpp:531
std::unordered_set< std::string > closed_this_step_
Definition: BlackoilWellModelGeneric.hpp:551
bool isOwner(const std::string &wname) const
Definition: BlackoilWellModelGeneric.hpp:274
void runWellPIScaling(const int reportStepIdx, DeferredLogger &local_deferredLogger)
std::vector< WellProdIndexCalculator< Scalar > > prod_index_calc_
Definition: BlackoilWellModelGeneric.hpp:547
Class for handling the blackoil well network model.
Definition: BlackoilWellModelNetworkGeneric.hpp:49
Class for handling the blackoil well model.
Definition: BlackoilWellModelWBP.hpp:42
Connection index mappings.
Definition: ConnectionIndexMap.hpp:33
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:53
void updateState(WellState< Scalar, IndexTraits > &well_state, GroupState< Scalar > &group_state)
Definition: GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: GasLiftGroupInfo.hpp:37
Definition: VFPProperties.hpp:40
Definition: WellInterfaceGeneric.hpp:53
Definition: WellState.hpp:66
ExcEnum
Definition: DeferredLogger.hpp:45
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilnewtonmethodparams.hpp:31
Definition: blackoilbioeffectsmodules.hh:43
Definition: BlackoilWellModelGeneric.hpp:87
static constexpr bool value
Definition: BlackoilWellModelGeneric.hpp:87
Definition: WGState.hpp:39