Schedule.hpp
Go to the documentation of this file.
1/*
2 Copyright 2013 Statoil ASA.
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef SCHEDULE_HPP
20#define SCHEDULE_HPP
21
22#include <map>
23#include <memory>
24
46
49
50
51/*
52 The DynamicState<std::shared_ptr<T>> pattern: The quantities in the Schedule
53 section like e.g. wellrates and completion properties are typically
54 characterized by the following behaviour:
55
56 1. They can be updated repeatedly at arbitrary points in the Schedule
57 section.
58
59 2. The value set at one timestep will apply until is explicitly set again at
60 a later timestep.
61
62 These properties are typically stored in a DynamicState<T> container; the
63 DynamicState<T> class is a container which implements this semantics:
64
65 1. It is legitimate to ask for an out-of-range value, you will then get the
66 last value which has been set.
67
68 2. When assigning an out-of-bounds value the container will append the
69 currently set value until correct length has been reached, and then the
70 new value will be assigned.
71
72 3. The DynamicState<T> has an awareness of the total length of the time
73 axis, trying to access values beyound that is illegal.
74
75 For many of the non-trival objects like eg Well and Group the DynamicState<>
76 contains a shared pointer to an underlying object, that way the fill operation
77 when the vector is resized is quite fast. The following pattern is quite
78 common for the Schedule implementation:
79
80
81 // Create a new well object.
82 std::shared_ptr<Well> new_well = this->getWell( well_name, time_step );
83
84 // Update the new well object with new settings from the deck, the
85 // updateXXXX() method will return true if the well object was actually
86 // updated:
87 if (new_well->updateRate( new_rate ))
88 this->dynamic_state.update( time_step, new_well);
89
90*/
91
92namespace Opm
93{
94 class Actions;
95 class Deck;
96 class DeckKeyword;
97 class DeckRecord;
98 class EclipseGrid;
99 class EclipseState;
100 class FieldPropsManager;
101 class Python;
102 class Runspec;
103 class RPTConfig;
104 class SCHEDULESection;
105 class SummaryState;
106 class TimeMap;
107 class UnitSystem;
108 class ErrorGuard;
109 class WListManager;
110 class UDQConfig;
111 class UDQActive;
112
113 class Schedule {
114 public:
117 using VFPProdMap = std::map<int, DynamicState<std::shared_ptr<VFPProdTable>>>;
118 using VFPInjMap = std::map<int, DynamicState<std::shared_ptr<VFPInjTable>>>;
119
120 Schedule() = default;
121 explicit Schedule(std::shared_ptr<const Python>) {}
122 Schedule(const Deck& deck,
123 const EclipseGrid& grid,
124 const FieldPropsManager& fp,
125 const Runspec &runspec,
126 const ParseContext& parseContext,
127 ErrorGuard& errors,
128 std::shared_ptr<const Python> python,
129 const RestartIO::RstState* rst = nullptr);
130
131 template<typename T>
132 Schedule(const Deck& deck,
133 const EclipseGrid& grid,
134 const FieldPropsManager& fp,
135 const Runspec &runspec,
136 const ParseContext& parseContext,
137 T&& errors,
138 std::shared_ptr<const Python> python,
139 const RestartIO::RstState* rst = nullptr);
140
141 Schedule(const Deck& deck,
142 const EclipseGrid& grid,
143 const FieldPropsManager& fp,
144 const Runspec &runspec,
145 std::shared_ptr<const Python> python,
146 const RestartIO::RstState* rst = nullptr);
147
148 Schedule(const Deck& deck,
149 const EclipseState& es,
150 const ParseContext& parseContext,
151 ErrorGuard& errors,
152 std::shared_ptr<const Python> python,
153 const RestartIO::RstState* rst = nullptr);
154
155 template <typename T>
156 Schedule(const Deck& deck,
157 const EclipseState& es,
158 const ParseContext& parseContext,
159 T&& errors,
160 std::shared_ptr<const Python> python,
161 const RestartIO::RstState* rst = nullptr);
162
163 Schedule(const Deck& deck,
164 const EclipseState& es,
165 std::shared_ptr<const Python> python,
166 const RestartIO::RstState* rst = nullptr);
167
168 // The constructor *without* the Python arg should really only be used from Python itself
169 Schedule(const Deck& deck,
170 const EclipseState& es,
171 const RestartIO::RstState* rst = nullptr);
172
174
175 /*
176 * If the input deck does not specify a start time, Eclipse's 1. Jan
177 * 1983 is defaulted
178 */
179 time_t getStartTime() const;
180 time_t posixStartTime() const;
181 time_t posixEndTime() const;
182 time_t simTime(size_t timeStep) const;
183 double seconds(size_t timeStep) const;
184 double stepLength(size_t timeStep) const;
185 std::optional<int> exitStatus() const;
186
187 const TimeMap& getTimeMap() const;
188
189 size_t numWells() const;
190 size_t numWells(size_t timestep) const;
191 bool hasWell(const std::string& wellName) const;
192 bool hasWell(const std::string& wellName, std::size_t timeStep) const;
193
194 std::vector<std::string> wellNames(const std::string& pattern, size_t timeStep, const std::vector<std::string>& matching_wells = {}) const;
195 std::vector<std::string> wellNames(const std::string& pattern) const;
196 std::vector<std::string> wellNames(size_t timeStep) const;
197 std::vector<std::string> wellNames() const;
198
199 std::vector<std::string> groupNames(const std::string& pattern, size_t timeStep) const;
200 std::vector<std::string> groupNames(size_t timeStep) const;
201 std::vector<std::string> groupNames(const std::string& pattern) const;
202 std::vector<std::string> groupNames() const;
203
204 void updateWell(std::shared_ptr<Well> well, size_t reportStep);
205 std::vector<std::string> changed_wells(size_t reportStep) const;
206 const Well& getWell(const std::string& wellName, size_t timeStep) const;
207 const Well& getWellatEnd(const std::string& well_name) const;
208 std::vector<Well> getWells(size_t timeStep) const;
209 std::vector<Well> getWellsatEnd() const;
210 void shut_well(const std::string& well_name, std::size_t report_step);
211 void stop_well(const std::string& well_name, std::size_t report_step);
212 void open_well(const std::string& well_name, std::size_t report_step);
213
214 std::vector<const Group*> getChildGroups2(const std::string& group_name, size_t timeStep) const;
215 std::vector<Well> getChildWells2(const std::string& group_name, size_t timeStep) const;
217 const Well::ProducerCMode& getGlobalWhistctlMmode(size_t timestep) const;
218
219 const UDQActive& udqActive(size_t timeStep) const;
220 const WellTestConfig& wtestConfig(size_t timestep) const;
221 const GConSale& gConSale(size_t timestep) const;
222 const GConSump& gConSump(size_t timestep) const;
223 const WListManager& getWListManager(size_t timeStep) const;
224 const UDQConfig& getUDQConfig(size_t timeStep) const;
225 const Action::Actions& actions(std::size_t timeStep) const;
226 void evalAction(const SummaryState& summary_state, size_t timeStep);
227
228 const RPTConfig& report_config(std::size_t timeStep) const;
229
230 GTNode groupTree(std::size_t report_step) const;
231 GTNode groupTree(const std::string& root_node, std::size_t report_step) const;
232 size_t numGroups() const;
233 size_t numGroups(size_t timeStep) const;
234 bool hasGroup(const std::string& groupName) const;
235 bool hasGroup(const std::string& groupName, std::size_t timeStep) const;
236 const Group& getGroup(const std::string& groupName, size_t timeStep) const;
237
238 const Tuning& getTuning(size_t timeStep) const;
241 void invalidNamePattern (const std::string& namePattern, std::size_t report_step, const ParseContext& parseContext, ErrorGuard& errors, const DeckKeyword& keyword) const;
242 const GuideRateConfig& guideRateConfig(size_t timeStep) const;
243
244 const RFTConfig& rftConfig() const;
245 const Events& getEvents() const;
246 const Events& getWellGroupEvents(const std::string& wellGroup) const;
247 bool hasWellGroupEvent(const std::string& wellGroup, uint64_t event_mask, size_t reportStep) const;
248 const Deck& getModifierDeck(size_t timeStep) const;
250 const VFPProdTable& getVFPProdTable(int table_id, size_t timeStep) const;
251 const VFPInjTable& getVFPInjTable(int table_id, size_t timeStep) const;
252 std::map<int, std::shared_ptr<const VFPProdTable> > getVFPProdTables(size_t timeStep) const;
253 std::map<int, std::shared_ptr<const VFPInjTable> > getVFPInjTables(size_t timeStep) const;
254 /*
255 Will remove all completions which are connected to cell which is not
256 active. Will scan through all wells and all timesteps.
257 */
259 size_t size() const;
260 const RestartConfig& restart() const;
262
263 void applyAction(size_t reportStep, const Action::ActionX& action, const Action::Result& result);
264 int getNupcol(size_t reportStep) const;
265
266 bool operator==(const Schedule& data) const;
267 std::shared_ptr<const Python> python() const;
268
269 /*
270 The cmp() function compares two schedule instances in a context aware
271 manner. Floating point numbers are compared with a tolerance. The
272 purpose of this comparison function is to implement regression tests
273 for the schedule instances created by loading a restart file.
274 */
275 static bool cmp(const Schedule& sched1, const Schedule& sched2, std::size_t report_step);
276
277 template<class Serializer>
278 void serializeOp(Serializer& serializer)
279 {
280 m_timeMap.serializeOp(serializer);
281 auto splitWells = splitDynMap(wells_static);
282 serializer.vector(splitWells.first);
283 serializer(splitWells.second);
284 auto splitGroups = splitDynMap(groups);
285 serializer.vector(splitGroups.first);
286 serializer(splitGroups.second);
287 m_oilvaporizationproperties.serializeOp(serializer);
288 m_events.serializeOp(serializer);
289 m_modifierDeck.serializeOp(serializer);
290 m_tuning.serializeOp(serializer);
291 m_messageLimits.serializeOp(serializer);
292 m_runspec.serializeOp(serializer);
293 auto splitvfpprod = splitDynMap<Map2>(vfpprod_tables);
294 serializer.vector(splitvfpprod.first);
295 serializer(splitvfpprod.second);
296 auto splitvfpinj = splitDynMap<Map2>(vfpinj_tables);
297 serializer.vector(splitvfpinj.first);
298 serializer(splitvfpinj.second);
299 wtest_config.serializeOp(serializer);
300 wlist_manager.serializeOp(serializer);
301 udq_config.serializeOp(serializer);
302 udq_active.serializeOp(serializer);
303 guide_rate_config.serializeOp(serializer);
304 gconsale.serializeOp(serializer);
305 gconsump.serializeOp(serializer);
306 global_whistctl_mode.template serializeOp<Serializer, false>(serializer);
307 m_actions.serializeOp(serializer);
308 rft_config.serializeOp(serializer);
309 m_nupcol.template serializeOp<Serializer, false>(serializer);
310 restart_config.serializeOp(serializer);
311 serializer.map(wellgroup_events);
312 if (!serializer.isSerializing()) {
313 reconstructDynMap(splitWells.first, splitWells.second, wells_static);
314 reconstructDynMap(splitGroups.first, splitGroups.second, groups);
315 reconstructDynMap<Map2>(splitvfpprod.first, splitvfpprod.second, vfpprod_tables);
316 reconstructDynMap<Map2>(splitvfpinj.first, splitvfpinj.second, vfpinj_tables);
317 }
318 }
319
320 private:
321 template<class Key, class Value> using Map2 = std::map<Key,Value>;
322 std::shared_ptr<const Python> python_handle;
323 TimeMap m_timeMap;
324 WellMap wells_static;
325 GroupMap groups;
326 DynamicState< OilVaporizationProperties > m_oilvaporizationproperties;
327 Events m_events;
328 DynamicVector< Deck > m_modifierDeck;
329 DynamicState<Tuning> m_tuning;
330 MessageLimits m_messageLimits;
331 Runspec m_runspec;
332 VFPProdMap vfpprod_tables;
333 VFPInjMap vfpinj_tables;
341 DynamicState<Well::ProducerCMode> global_whistctl_mode;
343 RFTConfig rft_config;
344 DynamicState<int> m_nupcol;
345 RestartConfig restart_config;
346 std::optional<int> exit_status;
347
348 std::map<std::string,Events> wellgroup_events;
349 void load_rst(const RestartIO::RstState& rst,
350 const EclipseGrid& grid,
351 const FieldPropsManager& fp,
352 const UnitSystem& unit_system);
353 void addWell(Well well, size_t report_step);
354 void addWell(const std::string& wellName,
355 const std::string& group,
356 int headI,
357 int headJ,
358 Phase preferredPhase,
359 double refDepth,
360 double drainageRadius,
361 bool allowCrossFlow,
362 bool automaticShutIn,
363 int pvt_table,
364 Well::GasInflowEquation gas_inflow,
365 size_t timeStep,
366 Connection::Order wellConnectionOrder,
367 const UnitSystem& unit_system);
368
370 GTNode groupTree(const std::string& root_node, std::size_t report_step, std::size_t level, const std::optional<std::string>& parent_name) const;
371 void updateGroup(std::shared_ptr<Group> group, size_t reportStep);
372 bool checkGroups(const ParseContext& parseContext, ErrorGuard& errors);
373 void updateUDQActive( std::size_t timeStep, std::shared_ptr<UDQActive> udq );
374 bool updateWellStatus( const std::string& well, size_t reportStep , Well::Status status, bool update_connections);
375 void addWellToGroup( const std::string& group_name, const std::string& well_name , size_t timeStep);
376 void iterateScheduleSection(std::shared_ptr<const Python> python, const std::string& input_path, const ParseContext& parseContext , ErrorGuard& errors, const SCHEDULESection& , const EclipseGrid& grid,
377 const FieldPropsManager& fp);
378 void addACTIONX(const Action::ActionX& action, std::size_t currentStep);
379 void addGroupToGroup( const std::string& parent_group, const std::string& child_group, size_t timeStep);
380 void addGroupToGroup( const std::string& parent_group, const Group& child_group, size_t timeStep);
381 void addGroup(const std::string& groupName , size_t timeStep, const UnitSystem& unit_system);
382 void addWell(const std::string& wellName, const DeckRecord& record, size_t timeStep, Connection::Order connection_order, const UnitSystem& unit_system);
383 void handleEXIT(const DeckKeyword& keyword , size_t report_step);
384 void handleUDQ(const DeckKeyword& keyword, size_t currentStep);
385 void handleWLIST(const DeckKeyword& keyword, size_t currentStep);
386 void handleCOMPORD(const ParseContext& parseContext, ErrorGuard& errors, const DeckKeyword& compordKeyword, size_t currentStep);
387 void handleWELSPECS( const SCHEDULESection&, size_t, size_t , const UnitSystem& unit_system, const ParseContext& parseContext, ErrorGuard& errors);
388 void handleWCONHIST( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
389 void handleWCONPROD( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
390 void handleWGRUPCON( const DeckKeyword& keyword, size_t currentStep);
391 void handleCOMPDAT( const DeckKeyword& keyword, size_t currentStep, const EclipseGrid& grid, const FieldPropsManager& fp, const ParseContext& parseContext, ErrorGuard& errors);
392 void handleCOMPLUMP( const DeckKeyword& keyword, size_t currentStep );
393 void handleWELSEGS( const DeckKeyword& keyword, size_t currentStep);
394 void handleCOMPSEGS( const DeckKeyword& keyword, size_t currentStep, const EclipseGrid& grid, const ParseContext& parseContext, ErrorGuard& errors);
395 void handleWSEGSICD( const DeckKeyword& keyword, size_t currentStep);
396 // TODO: we should incorporate ParseContext and ErrorGuard, including the above keyword
397 void handleWSEGVALV( const DeckKeyword& keyword, size_t currentStep);
398
399 void handleWCONINJE( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
400 void handleWFOAM( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
401 void handleWPOLYMER( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
402 void handleWSALT( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
403 void handleWSOLVENT( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
404 void handleWTRACER( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
405 void handleWTEMP( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
406 void handleWPMITAB( const DeckKeyword& keyword, const size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
407 void handleWSKPTAB( const DeckKeyword& keyword, const size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
408 void handleWINJTEMP( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
409 void handleWCONINJH(const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
410 void handleWELOPEN( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors, const std::vector<std::string>& matching_wells = {});
411 void handleWELTARG( const SCHEDULESection&, const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
412 void handleGCONINJE( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
413 void handleGCONPROD( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
414 void handleGEFAC( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
415 void handleGCONSALE( const DeckKeyword& keyword, size_t currentStep, const UnitSystem& unit_system);
416 void handleGCONSUMP( const DeckKeyword& keyword, size_t currentStep, const UnitSystem& unit_system);
417 void handleGUIDERAT( const DeckKeyword& keyword, size_t currentStep);
418 void handleLINCOM( const DeckKeyword& keyword, size_t currentStep);
419 void handleWEFAC( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
420
421 void handleTUNING( const DeckKeyword& keyword, size_t currentStep);
422 void handlePYACTION( std::shared_ptr<const Python> python, const std::string& input_path, const DeckKeyword& keyword, size_t currentStep);
423 void handleNUPCOL( const DeckKeyword& keyword, size_t currentStep);
424 void handleGRUPTREE( const DeckKeyword& keyword, size_t currentStep, const UnitSystem& unit_system, const ParseContext& parseContext, ErrorGuard& errors);
425 void handleGRUPNET( const DeckKeyword& keyword, size_t currentStep, const UnitSystem& unit_system);
426 void handleWRFT( const DeckKeyword& keyword, size_t currentStep);
427 void handleWTEST( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
428 void handleWRFTPLT( const DeckKeyword& keyword, size_t currentStep);
429 void handleWPIMULT( const DeckKeyword& keyword, size_t currentStep);
430 void handleDRSDT( const DeckKeyword& keyword, size_t currentStep);
431 void handleDRVDT( const DeckKeyword& keyword, size_t currentStep);
432 void handleDRSDTR( const DeckKeyword& keyword, size_t currentStep);
433 void handleDRVDTR( const DeckKeyword& keyword, size_t currentStep);
434 void handleVAPPARS( const DeckKeyword& keyword, size_t currentStep);
435 void handleWECON( const DeckKeyword& keyword, size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
436 void handleWHISTCTL(const DeckKeyword& keyword, std::size_t currentStep, const ParseContext& parseContext, ErrorGuard& errors);
437 void handleMESSAGES(const DeckKeyword& keyword, size_t currentStep);
438 void handleRPTSCHED(const DeckKeyword& keyword, size_t currentStep);
439 void handleVFPPROD(const DeckKeyword& vfpprodKeyword, const UnitSystem& unit_system, size_t currentStep);
440 void handleVFPINJ(const DeckKeyword& vfpprodKeyword, const UnitSystem& unit_system, size_t currentStep);
441 void checkUnhandledKeywords( const SCHEDULESection& ) const;
442 void checkIfAllConnectionsIsShut(size_t currentStep);
443 void handleKeyword(std::shared_ptr<const Python> python,
444 const std::string& input_path,
445 size_t currentStep,
446 const SCHEDULESection& section,
447 size_t keywordIdx,
448 const DeckKeyword& keyword,
449 const ParseContext& parseContext, ErrorGuard& errors,
450 const EclipseGrid& grid,
451 const FieldPropsManager& fp,
452 const UnitSystem& unit_system,
453 std::vector<std::pair<const DeckKeyword*, size_t > >& rftProperties);
454 void addWellGroupEvent(const std::string& wellGroup, ScheduleEvents::Events event, size_t reportStep);
455
456 template<template<class, class> class Map, class Type, class Key>
457 std::pair<std::vector<Type>, std::vector<std::pair<Key, std::vector<size_t>>>>
458 splitDynMap(const Map<Key, Opm::DynamicState<Type>>& map)
459 {
460 // we have to pack the unique ptrs separately, and use an index map
461 // to allow reconstructing the appropriate structures.
462 std::vector<std::pair<Key, std::vector<size_t>>> asMap;
463 std::vector<Type> unique;
464 for (const auto& it : map) {
465 auto indices = it.second.split(unique);
466 asMap.push_back(std::make_pair(it.first, indices));
467 }
468
469 return std::make_pair(unique, asMap);
470 }
471
472 template<template<class, class> class Map, class Type, class Key>
473 void reconstructDynMap(const std::vector<Type>& unique,
474 const std::vector<std::pair<Key, std::vector<size_t>>>& asMap,
475 Map<Key, Opm::DynamicState<Type>>& result)
476 {
477 for (const auto& it : asMap) {
478 result[it.first].reconstruct(unique, it.second);
479 }
480 }
481 };
482}
483
484#endif
const char *const string
Definition: cJSON.h:170
Definition: ActionX.hpp:64
Definition: Actions.hpp:39
Definition: ActionResult.hpp:89
Simple class capturing active cells of a grid.
Definition: ActiveGridCells.hpp:35
Order
Definition: parser/eclipse/EclipseState/Schedule/Well/connection.hpp:67
Definition: DeckKeyword.hpp:38
Definition: Deck.hpp:115
Definition: DeckRecord.hpp:32
Definition: DynamicState.hpp:58
Definition: DynamicVector.hpp:38
Definition: EclipseGrid.hpp:54
Definition: EclipseState.hpp:63
Definition: ErrorGuard.hpp:29
Definition: Events.hpp:122
void serializeOp(Serializer &serializer)
Definition: Events.hpp:135
Definition: FieldPropsManager.hpp:33
Definition: GConSale.hpp:32
Definition: GConSump.hpp:32
Definition: GTNode.hpp:31
Definition: parser/eclipse/EclipseState/Schedule/Group/group.hpp:36
Definition: GuideRateConfig.hpp:34
Definition: MessageLimits.hpp:101
void serializeOp(Serializer &serializer)
Definition: MessageLimits.hpp:146
Definition: OilVaporizationProperties.hpp:34
Definition: ParseContext.hpp:84
Definition: RFTConfig.hpp:33
void serializeOp(Serializer &serializer)
Definition: RFTConfig.hpp:83
Definition: RPTConfig.hpp:29
Definition: RestartConfig.hpp:326
void serializeOp(Serializer &serializer)
Definition: RestartConfig.hpp:354
Definition: Runspec.hpp:264
void serializeOp(Serializer &serializer)
Definition: Runspec.hpp:285
Definition: DeckSection.hpp:115
Definition: Schedule.hpp:113
Schedule(const Deck &deck, const EclipseGrid &grid, const FieldPropsManager &fp, const Runspec &runspec, std::shared_ptr< const Python > python, const RestartIO::RstState *rst=nullptr)
const TimeMap & getTimeMap() const
std::vector< std::string > changed_wells(size_t reportStep) const
void filterConnections(const ActiveGridCells &grid)
void serializeOp(Serializer &serializer)
Definition: Schedule.hpp:278
const GConSale & gConSale(size_t timestep) const
const RPTConfig & report_config(std::size_t timeStep) const
std::vector< std::string > wellNames(const std::string &pattern, size_t timeStep, const std::vector< std::string > &matching_wells={}) const
Schedule()=default
std::shared_ptr< const Python > python() const
const DynamicState< Tuning > & getTuning() const
bool hasWell(const std::string &wellName, std::size_t timeStep) const
const VFPProdTable & getVFPProdTable(int table_id, size_t timeStep) const
const WellTestConfig & wtestConfig(size_t timestep) const
std::map< int, DynamicState< std::shared_ptr< VFPInjTable > > > VFPInjMap
Definition: Schedule.hpp:118
std::vector< Well > getWellsatEnd() const
void applyAction(size_t reportStep, const Action::ActionX &action, const Action::Result &result)
const MessageLimits & getMessageLimits() const
bool hasGroup(const std::string &groupName, std::size_t timeStep) const
const UDQConfig & getUDQConfig(size_t timeStep) const
void updateWell(std::shared_ptr< Well > well, size_t reportStep)
void open_well(const std::string &well_name, std::size_t report_step)
const Events & getWellGroupEvents(const std::string &wellGroup) const
RestartConfig & restart()
const VFPInjTable & getVFPInjTable(int table_id, size_t timeStep) const
std::map< int, std::shared_ptr< const VFPInjTable > > getVFPInjTables(size_t timeStep) const
std::vector< std::string > wellNames() const
size_t numWells() const
Schedule(const Deck &deck, const EclipseGrid &grid, const FieldPropsManager &fp, const Runspec &runspec, const ParseContext &parseContext, ErrorGuard &errors, std::shared_ptr< const Python > python, const RestartIO::RstState *rst=nullptr)
int getNupcol(size_t reportStep) const
Schedule(const Deck &deck, const EclipseState &es, std::shared_ptr< const Python > python, const RestartIO::RstState *rst=nullptr)
const GuideRateConfig & guideRateConfig(size_t timeStep) const
GTNode groupTree(std::size_t report_step) const
void invalidNamePattern(const std::string &namePattern, std::size_t report_step, const ParseContext &parseContext, ErrorGuard &errors, const DeckKeyword &keyword) const
time_t simTime(size_t timeStep) const
const GConSump & gConSump(size_t timestep) const
time_t posixEndTime() const
bool hasGroup(const std::string &groupName) const
bool hasOilVaporizationProperties() const
std::vector< std::string > wellNames(size_t timeStep) const
const WListManager & getWListManager(size_t timeStep) const
Schedule(const Deck &deck, const EclipseState &es, const ParseContext &parseContext, T &&errors, std::shared_ptr< const Python > python, const RestartIO::RstState *rst=nullptr)
const RFTConfig & rftConfig() const
const Group & getGroup(const std::string &groupName, size_t timeStep) const
std::vector< std::string > wellNames(const std::string &pattern) const
Schedule(const Deck &deck, const EclipseState &es, const RestartIO::RstState *rst=nullptr)
Schedule(const Deck &deck, const EclipseState &es, const ParseContext &parseContext, ErrorGuard &errors, std::shared_ptr< const Python > python, const RestartIO::RstState *rst=nullptr)
void evalAction(const SummaryState &summary_state, size_t timeStep)
const UDQActive & udqActive(size_t timeStep) const
size_t numGroups() const
bool operator==(const Schedule &data) const
const Well::ProducerCMode & getGlobalWhistctlMmode(size_t timestep) const
void shut_well(const std::string &well_name, std::size_t report_step)
void stop_well(const std::string &well_name, std::size_t report_step)
const Events & getEvents() const
size_t size() const
GTNode groupTree(const std::string &root_node, std::size_t report_step) const
OrderedMap< std::string, DynamicState< std::shared_ptr< Well > > > WellMap
Definition: Schedule.hpp:115
const RestartConfig & restart() const
static bool cmp(const Schedule &sched1, const Schedule &sched2, std::size_t report_step)
std::vector< std::string > groupNames(size_t timeStep) const
const Tuning & getTuning(size_t timeStep) const
std::map< int, std::shared_ptr< const VFPProdTable > > getVFPProdTables(size_t timeStep) const
size_t numGroups(size_t timeStep) const
Schedule(std::shared_ptr< const Python >)
Definition: Schedule.hpp:121
time_t getStartTime() const
std::optional< int > exitStatus() const
std::vector< Well > getWells(size_t timeStep) const
double seconds(size_t timeStep) const
double stepLength(size_t timeStep) const
const Action::Actions & actions(std::size_t timeStep) const
bool hasWell(const std::string &wellName) const
size_t numWells(size_t timestep) const
const Well & getWellatEnd(const std::string &well_name) const
OrderedMap< std::string, DynamicState< std::shared_ptr< Group > > > GroupMap
Definition: Schedule.hpp:116
std::vector< Well > getChildWells2(const std::string &group_name, size_t timeStep) const
std::vector< const Group * > getChildGroups2(const std::string &group_name, size_t timeStep) const
std::map< int, DynamicState< std::shared_ptr< VFPProdTable > > > VFPProdMap
Definition: Schedule.hpp:117
std::vector< std::string > groupNames(const std::string &pattern) const
Schedule(const Deck &deck, const EclipseGrid &grid, const FieldPropsManager &fp, const Runspec &runspec, const ParseContext &parseContext, T &&errors, std::shared_ptr< const Python > python, const RestartIO::RstState *rst=nullptr)
time_t posixStartTime() const
std::vector< std::string > groupNames(const std::string &pattern, size_t timeStep) const
std::vector< std::string > groupNames() const
const Well & getWell(const std::string &wellName, size_t timeStep) const
bool hasWellGroupEvent(const std::string &wellGroup, uint64_t event_mask, size_t reportStep) const
static Schedule serializeObject()
const OilVaporizationProperties & getOilVaporizationProperties(size_t timestep) const
const Deck & getModifierDeck(size_t timeStep) const
Definition: Serializer.hpp:38
Definition: SummaryState.hpp:65
Definition: TimeMap.hpp:40
void serializeOp(Serializer &serializer)
Definition: TimeMap.hpp:81
Definition: UDQActive.hpp:35
Definition: UDQConfig.hpp:41
Definition: UnitSystem.hpp:32
Definition: VFPInjTable.hpp:33
Definition: VFPProdTable.hpp:36
Definition: WListManager.hpp:28
Definition: custom-opm-common/opm-common/opm/parser/eclipse/EclipseState/Schedule/Well/well.hpp:58
Status
Definition: custom-opm-common/opm-common/opm/parser/eclipse/EclipseState/Schedule/Well/well.hpp:61
ProducerCMode
Definition: custom-opm-common/opm-common/opm/parser/eclipse/EclipseState/Schedule/Well/well.hpp:101
GasInflowEquation
Definition: custom-opm-common/opm-common/opm/parser/eclipse/EclipseState/Schedule/Well/well.hpp:152
Definition: WellTestConfig.hpp:29
Events
Definition: Events.hpp:30
UDAKeyword keyword(UDAControl control)
std::vector< typename std::result_of< F(typename C::const_iterator::value_type &) >::type > map(F f, const C &src)
Definition: Functional.hpp:84
status
Definition: value_status.hpp:27
Definition: A.hpp:4
Phase
Definition: Runspec.hpp:34
static std::string data()
Definition: exprtk.hpp:40022
Definition: state.hpp:37
Definition: Tuning.hpp:27