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