My Project
Loading...
Searching...
No Matches
Schedule.hpp
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 <cstddef>
23#include <ctime>
24#include <functional>
25#include <iosfwd>
26#include <map>
27#include <memory>
28#include <optional>
29#include <set>
30#include <string>
31#include <unordered_map>
32#include <utility>
33#include <vector>
34
35#include <opm/input/eclipse/Schedule/Action/SimulatorUpdate.hpp>
36#include <opm/input/eclipse/Schedule/Action/WGNames.hpp>
37#include <opm/input/eclipse/Schedule/CompletedCells.hpp>
38#include <opm/input/eclipse/Schedule/Group/Group.hpp>
39#include <opm/input/eclipse/Schedule/ScheduleDeck.hpp>
40#include <opm/input/eclipse/Schedule/ScheduleState.hpp>
41#include <opm/input/eclipse/Schedule/ScheduleStatic.hpp>
42#include <opm/input/eclipse/Schedule/Well/PAvg.hpp>
43#include <opm/input/eclipse/Schedule/Well/Well.hpp>
44#include <opm/input/eclipse/Schedule/WriteRestartFileEvents.hpp>
45#include <opm/input/eclipse/Units/UnitSystem.hpp>
46
47namespace Opm
48{
49 namespace Action {
50 class ActionX;
51 class PyAction;
52 class State;
53 }
54 class ActiveGridCells;
55 class Deck;
56 class DeckKeyword;
57 class DeckRecord;
58 enum class ConnectionOrder;
59 class EclipseGrid;
60 class EclipseState;
61 class ErrorGuard;
62 class FieldPropsManager;
63 class GasLiftOpt;
64 class GTNode;
65 class GuideRateConfig;
66 class GuideRateModel;
67 class HandlerContext;
68 enum class InputErrorAction;
69 class ParseContext;
70 class Python;
71 class Runspec;
72 class RPTConfig;
73 class ScheduleGrid;
74 class SCHEDULESection;
75 class SegmentMatcher;
76 class SummaryState;
77 class TracerConfig;
78 class UDQConfig;
79 class Well;
80 enum class WellGasInflowEquation;
81 class WellMatcher;
82 enum class WellProducerCMode;
83 enum class WellStatus;
84 class WelSegsSet;
85 class WellTestConfig;
86
87 namespace RestartIO { struct RstState; }
88
89 class Schedule {
90 public:
91 Schedule() = default;
92
93 explicit Schedule(std::shared_ptr<const Python> python_handle);
94
108 Schedule(const Deck& deck,
109 const EclipseGrid& grid,
110 const FieldPropsManager& fp,
111 const Runspec &runspec,
112 const ParseContext& parseContext,
113 ErrorGuard& errors,
114 std::shared_ptr<const Python> python,
115 const bool lowActionParsingStrictness = false,
116 const bool keepKeywords = true,
117 const std::optional<int>& output_interval = {},
118 const RestartIO::RstState* rst = nullptr,
119 const TracerConfig* tracer_config = nullptr);
120
121 template<typename T>
122 Schedule(const Deck& deck,
123 const EclipseGrid& grid,
124 const FieldPropsManager& fp,
125 const Runspec &runspec,
126 const ParseContext& parseContext,
127 T&& errors,
128 std::shared_ptr<const Python> python,
129 const bool lowActionParsingStrictness = false,
130 const bool keepKeywords = true,
131 const std::optional<int>& output_interval = {},
132 const RestartIO::RstState* rst = nullptr,
133 const TracerConfig* tracer_config = nullptr);
134
135 Schedule(const Deck& deck,
136 const EclipseGrid& grid,
137 const FieldPropsManager& fp,
138 const Runspec &runspec,
139 std::shared_ptr<const Python> python,
140 const bool lowActionParsingStrictness = false,
141 const bool keepKeywords = true,
142 const std::optional<int>& output_interval = {},
143 const RestartIO::RstState* rst = nullptr,
144 const TracerConfig* tracer_config = nullptr);
145
146 Schedule(const Deck& deck,
147 const EclipseState& es,
148 const ParseContext& parseContext,
149 ErrorGuard& errors,
150 std::shared_ptr<const Python> python,
151 const bool lowActionParsingStrictness = false,
152 const bool keepKeywords = true,
153 const std::optional<int>& output_interval = {},
154 const RestartIO::RstState* rst = nullptr);
155
156 template <typename T>
157 Schedule(const Deck& deck,
158 const EclipseState& es,
159 const ParseContext& parseContext,
160 T&& errors,
161 std::shared_ptr<const Python> python,
162 const bool lowActionParsingStrictness = false,
163 const bool keepKeywords = true,
164 const std::optional<int>& output_interval = {},
165 const RestartIO::RstState* rst = nullptr);
166
167 Schedule(const Deck& deck,
168 const EclipseState& es,
169 std::shared_ptr<const Python> python,
170 const bool lowActionParsingStrictness = false,
171 const bool keepKeywords = true,
172 const std::optional<int>& output_interval = {},
173 const RestartIO::RstState* rst = nullptr);
174
175 // The constructor *without* the Python arg should really only be used from Python itself
176 Schedule(const Deck& deck,
177 const EclipseState& es,
178 const std::optional<int>& output_interval = {},
179 const RestartIO::RstState* rst = nullptr);
180
181 ~Schedule() = default;
182
183 static Schedule serializationTestObject();
184
185 /*
186 * If the input deck does not specify a start time, Eclipse's 1. Jan
187 * 1983 is defaulted
188 */
189 std::time_t getStartTime() const;
190 std::time_t posixStartTime() const;
191 std::time_t posixEndTime() const;
192 std::time_t simTime(std::size_t timeStep) const;
193 double seconds(std::size_t timeStep) const;
194 double stepLength(std::size_t timeStep) const;
195 std::optional<int> exitStatus() const;
196 const UnitSystem& getUnits() const { return this->m_static.m_unit_system; }
197 const Runspec& runspec() const { return this->m_static.m_runspec; }
198
199 std::size_t numWells() const;
200 std::size_t numWells(std::size_t timestep) const;
201 bool hasWell(const std::string& wellName) const;
202 bool hasWell(const std::string& wellName, std::size_t timeStep) const;
203
204 WellMatcher wellMatcher(std::size_t report_step) const;
205 std::function<std::unique_ptr<SegmentMatcher>()> segmentMatcherFactory(std::size_t report_step) const;
206 std::vector<std::string> wellNames(const std::string& pattern, std::size_t timeStep, const std::vector<std::string>& matching_wells = {}) const;
207 std::vector<std::string> wellNames(const std::string& pattern) const;
208 std::vector<std::string> wellNames(std::size_t timeStep) const;
209 std::vector<std::string> wellNames() const;
210
211 bool hasGroup(const std::string& groupName, std::size_t timeStep) const;
212 std::vector<std::string> groupNames(const std::string& pattern, std::size_t timeStep) const;
213 std::vector<std::string> groupNames(std::size_t timeStep) const;
214 std::vector<std::string> groupNames(const std::string& pattern) const;
215 std::vector<std::string> groupNames() const;
216 /*
217 The restart_groups function returns a vector of groups pointers which
218 is organized as follows:
219
220 1. The number of elements is WELLDIMS::MAXGROUPS + 1
221 2. The elements are sorted according to group.insert_index().
222 3. If there are less than WELLDIMS::MAXGROUPS nullptr is used.
223 4. The very last element corresponds to the FIELD group.
224 */
225 std::vector<const Group*> restart_groups(std::size_t timeStep) const;
226
227 std::vector<std::string> changed_wells(std::size_t reportStep) const;
228 const Well& getWell(std::size_t well_index, std::size_t timeStep) const;
229 const Well& getWell(const std::string& wellName, std::size_t timeStep) const;
230 const Well& getWellatEnd(const std::string& well_name) const;
231 // get the list of the constant flux aquifer specified in the whole schedule
232 std::unordered_set<int> getAquiferFluxSchedule() const;
233 std::vector<Well> getWells(std::size_t timeStep) const;
234 std::vector<Well> getWellsatEnd() const;
235
236 const std::unordered_map<std::string, std::set<int>>& getPossibleFutureConnections() const;
237
238 void shut_well(const std::string& well_name, std::size_t report_step);
239 void shut_well(const std::string& well_name);
240 void stop_well(const std::string& well_name, std::size_t report_step);
241 void stop_well(const std::string& well_name);
242 void open_well(const std::string& well_name, std::size_t report_step);
243 void open_well(const std::string& well_name);
244 void clear_event(ScheduleEvents::Events, std::size_t report_step);
245 void applyWellProdIndexScaling(const std::string& well_name, const std::size_t reportStep, const double scalingFactor);
246
247 std::vector<const Group*> getChildGroups2(const std::string& group_name, std::size_t timeStep) const;
248 std::vector<Well> getChildWells2(const std::string& group_name, std::size_t timeStep) const;
249 WellProducerCMode getGlobalWhistctlMmode(std::size_t timestep) const;
250
251 const UDQConfig& getUDQConfig(std::size_t timeStep) const;
252 void evalAction(const SummaryState& summary_state, std::size_t timeStep);
253
254 GTNode groupTree(std::size_t report_step) const;
255 GTNode groupTree(const std::string& root_node, std::size_t report_step) const;
256 const Group& getGroup(const std::string& groupName, std::size_t timeStep) const;
257
258 std::optional<std::size_t> first_RFT() const;
259 /*
260 Will remove all completions which are connected to cell which is not
261 active. Will scan through all wells and all timesteps.
262 */
263 void filterConnections(const ActiveGridCells& grid);
264 std::size_t size() const;
265
266 bool write_rst_file(std::size_t report_step) const;
267 const std::map< std::string, int >& rst_keywords( size_t timestep ) const;
268
269 /*
270 The applyAction() is invoked from the simulator *after* an ACTIONX has
271 evaluated to true. The return value is a small structure with
272 'information' which the simulator should take into account when
273 updating internal datastructures after the ACTIONX keywords have been
274 applied.
275 */
276 SimulatorUpdate applyAction(std::size_t reportStep,
277 const Action::ActionX& action,
278 const std::vector<std::string>& matching_wells,
279 const std::unordered_map<std::string, double>& wellpi);
280
281 SimulatorUpdate applyAction(std::size_t reportStep,
282 const Action::ActionX& action,
283 const std::vector<std::string>& matching_wells,
284 const std::unordered_map<std::string, float>& wellpi);
285 /*
286 The runPyAction() will run the Python script in a PYACTION keyword. In
287 the case of Schedule updates the recommended way of doing that from
288 PYACTION is to invoke a "normal" ACTIONX keyword internally from the
289 Python code. he return value from runPyAction() comes from such a
290 internal ACTIONX.
291 */
292 SimulatorUpdate runPyAction(std::size_t reportStep, const Action::PyAction& pyaction, Action::State& action_state, EclipseState& ecl_state, SummaryState& summary_state);
293
294
295 const GasLiftOpt& glo(std::size_t report_step) const;
296
297 bool operator==(const Schedule& data) const;
298 std::shared_ptr<const Python> python() const;
299
300
301 const ScheduleState& back() const;
302 const ScheduleState& operator[](std::size_t index) const;
303 std::vector<ScheduleState>::const_iterator begin() const;
304 std::vector<ScheduleState>::const_iterator end() const;
305 void create_next(const time_point& start_time, const std::optional<time_point>& end_time);
306 void create_next(const ScheduleBlock& block);
307 void create_first(const time_point& start_time, const std::optional<time_point>& end_time);
308
309 void treat_critical_as_non_critical(bool value) { this->m_treat_critical_as_non_critical = value; }
310
311 /*
312 The cmp() function compares two schedule instances in a context aware
313 manner. Floating point numbers are compared with a tolerance. The
314 purpose of this comparison function is to implement regression tests
315 for the schedule instances created by loading a restart file.
316 */
317 static bool cmp(const Schedule& sched1, const Schedule& sched2, std::size_t report_step);
318 void applyKeywords(std::vector<std::unique_ptr<DeckKeyword>>& keywords, std::size_t report_step);
319 void applyKeywords(std::vector<std::unique_ptr<DeckKeyword>>& keywords);
320
321 template<class Serializer>
322 void serializeOp(Serializer& serializer)
323 {
324 serializer(this->m_static);
325 serializer(this->m_sched_deck);
326 serializer(this->action_wgnames);
327 serializer(this->exit_status);
328 serializer(this->snapshots);
329 serializer(this->restart_output);
330 serializer(this->completed_cells);
331 serializer(this->m_treat_critical_as_non_critical);
332 serializer(this->current_report_step);
333 serializer(this->m_lowActionParsingStrictness);
334 serializer(this->simUpdateFromPython);
335
336 this->template pack_unpack<PAvg>(serializer);
337 this->template pack_unpack<WellTestConfig>(serializer);
338 this->template pack_unpack<GConSale>(serializer);
339 this->template pack_unpack<GConSump>(serializer);
340 this->template pack_unpack<GroupEconProductionLimits>(serializer);
341 this->template pack_unpack<WListManager>(serializer);
342 this->template pack_unpack<Network::ExtNetwork>(serializer);
343 this->template pack_unpack<Network::Balance>(serializer);
344 this->template pack_unpack<RPTConfig>(serializer);
345 this->template pack_unpack<Action::Actions>(serializer);
346 this->template pack_unpack<UDQActive>(serializer);
347 this->template pack_unpack<UDQConfig>(serializer);
348 this->template pack_unpack<NameOrder>(serializer);
349 this->template pack_unpack<GroupOrder>(serializer);
350 this->template pack_unpack<GuideRateConfig>(serializer);
351 this->template pack_unpack<GasLiftOpt>(serializer);
352 this->template pack_unpack<RFTConfig>(serializer);
353 this->template pack_unpack<RSTConfig>(serializer);
354 this->template pack_unpack<ScheduleState::BHPDefaults>(serializer);
355 this->template pack_unpack<Source>(serializer);
356
357 this->template pack_unpack_map<int, VFPProdTable>(serializer);
358 this->template pack_unpack_map<int, VFPInjTable>(serializer);
359 this->template pack_unpack_map<std::string, Group>(serializer);
360 this->template pack_unpack_map<std::string, Well>(serializer);
361
362 // If we are deserializing we need to setup the pointer to the
363 // unit system since this is process specific. This is safe
364 // because we set the same value in all well instances.
365 // We do some redundant assignments as these are shared_ptr's
366 // with multiple pointers to any given instance, but it is not
367 // significant so let's keep it simple.
368 if (!serializer.isSerializing()) {
369 for (auto& snapshot : snapshots) {
370 for (auto& well : snapshot.wells) {
371 well.second->updateUnitSystem(&m_static.m_unit_system);
372 }
373 }
374 }
375 }
376
377 template <typename T, class Serializer>
378 void pack_unpack(Serializer& serializer) {
379 std::vector<T> value_list;
380 std::vector<std::size_t> index_list;
381
382 if (serializer.isSerializing())
383 this->template pack_state<T>(value_list, index_list);
384
385 serializer(value_list);
386 serializer(index_list);
387
388 if (!serializer.isSerializing())
389 this->template unpack_state<T>(value_list, index_list);
390 }
391
392 template <typename T>
393 std::vector<std::pair<std::size_t, T>> unique() const {
394 std::vector<std::pair<std::size_t, T>> values;
395 for (std::size_t index = 0; index < this->snapshots.size(); index++) {
396 const auto& member = this->snapshots[index].get<T>();
397 const auto& value = member.get();
398 if (values.empty() || !(value == values.back().second))
399 values.push_back( std::make_pair(index, value));
400 }
401 return values;
402 }
403
404
405 template <typename T>
406 void pack_state(std::vector<T>& value_list, std::vector<std::size_t>& index_list) const {
407 auto unique_values = this->template unique<T>();
408 for (auto& [index, value] : unique_values) {
409 value_list.push_back( std::move(value) );
410 index_list.push_back( index );
411 }
412 }
413
414
415 template <typename T>
416 void unpack_state(const std::vector<T>& value_list, const std::vector<std::size_t>& index_list) {
417 std::size_t unique_index = 0;
418 while (unique_index < value_list.size()) {
419 const auto& value = value_list[unique_index];
420 const auto& first_index = index_list[unique_index];
421 auto last_index = this->snapshots.size();
422 if (unique_index < (value_list.size() - 1))
423 last_index = index_list[unique_index + 1];
424
425 auto& target_state = this->snapshots[first_index];
426 target_state.get<T>().update( std::move(value) );
427 for (std::size_t index=first_index + 1; index < last_index; index++)
428 this->snapshots[index].get<T>().update( target_state.get<T>() );
429
430 unique_index++;
431 }
432 }
433
434
435 template <typename K, typename T, class Serializer>
436 void pack_unpack_map(Serializer& serializer) {
437 std::vector<T> value_list;
438 std::vector<std::size_t> index_list;
439
440 if (serializer.isSerializing())
441 pack_map<K,T>(value_list, index_list);
442
443 serializer(value_list);
444 serializer(index_list);
445
446 if (!serializer.isSerializing())
447 unpack_map<K,T>(value_list, index_list);
448 }
449
450
451 template <typename K, typename T>
452 void pack_map(std::vector<T>& value_list,
453 std::vector<std::size_t>& index_list) {
454
455 const auto& last_map = this->snapshots.back().get_map<K,T>();
456 std::vector<K> key_list{ last_map.keys() };
457 std::unordered_map<K,T> current_value;
458
459 for (std::size_t index = 0; index < this->snapshots.size(); index++) {
460 auto& state = this->snapshots[index];
461 const auto& current_map = state.template get_map<K,T>();
462 for (const auto& key : key_list) {
463 auto& value = current_map.get_ptr(key);
464 if (value) {
465 auto it = current_value.find(key);
466 if (it == current_value.end() || !(*value == it->second)) {
467 value_list.push_back( *value );
468 index_list.push_back( index );
469
470 current_value[key] = *value;
471 }
472 }
473 }
474 }
475 }
476
477
478 template <typename K, typename T>
479 void unpack_map(const std::vector<T>& value_list,
480 const std::vector<std::size_t>& index_list) {
481
482 std::unordered_map<K, std::vector<std::pair<std::size_t, T>>> storage;
483 for (std::size_t storage_index = 0; storage_index < value_list.size(); storage_index++) {
484 const auto& value = value_list[storage_index];
485 const auto& time_index = index_list[storage_index];
486
487 storage[ value.name() ].emplace_back( time_index, value );
488 }
489
490 for (const auto& [key, values] : storage) {
491 for (std::size_t unique_index = 0; unique_index < values.size(); unique_index++) {
492 const auto& [time_index, value] = values[unique_index];
493 auto last_index = this->snapshots.size();
494 if (unique_index < (values.size() - 1))
495 last_index = values[unique_index + 1].first;
496
497 auto& map_value = this->snapshots[time_index].template get_map<K,T>();
498 map_value.update(std::move(value));
499
500 for (std::size_t index=time_index + 1; index < last_index; index++) {
501 auto& forward_map = this->snapshots[index].template get_map<K,T>();
502 forward_map.update( key, map_value );
503 }
504 }
505 }
506 }
507
508 friend std::ostream& operator<<(std::ostream& os, const Schedule& sched);
509 void dump_deck(std::ostream& os) const;
510
511 private:
512 friend class HandlerContext;
513
514 // Please update the member functions
515 // - operator==(const Schedule&) const
516 // - serializationTestObject()
517 // - serializeOp(Serializer&)
518 // when you update/change this list of data members.
519 bool m_treat_critical_as_non_critical = false;
520 ScheduleStatic m_static{};
521 ScheduleDeck m_sched_deck{};
522 Action::WGNames action_wgnames{};
523 std::optional<int> exit_status{};
524 std::vector<ScheduleState> snapshots{};
525 WriteRestartFileEvents restart_output{};
526 CompletedCells completed_cells{};
527
528 // Boolean indicating the strictness of parsing process for ActionX and PyAction.
529 // If lowActionParsingStrictness is true, the simulator tries to apply unsupported
530 // keywords, if lowActionParsingStrictness is false, the simulator only applies
531 // supported keywords.
532 bool m_lowActionParsingStrictness = false;
533
534 // This unordered_map contains possible future connections of wells that might get added through an ACTIONX.
535 // For parallel runs, this unordered_map is retrieved by the grid partitioner to ensure these connections
536 // end up on the same partition.
537 std::unordered_map<std::string, std::set<int>> possibleFutureConnections;
538
539 // The current_report_step is set to the current report step when a PYACTION call is executed.
540 // This is needed since the Schedule object does not know the current report step of the simulator and
541 // we only allow PYACTIONS for the current and future report steps.
542 std::size_t current_report_step = 0;
543 // The simUpdateFromPython points to a SimulatorUpdate collecting all updates from one PYACTION call.
544 // The SimulatorUpdate is reset before a new PYACTION call is executed.
545 // It is a shared_ptr, so a Schedule can be constructed using the copy constructor sharing the simUpdateFromPython.
546 // The copy constructor is needed for creating a mocked simulator (msim).
547 std::shared_ptr<SimulatorUpdate> simUpdateFromPython{};
548
549 void load_rst(const RestartIO::RstState& rst,
550 const TracerConfig& tracer_config,
551 const ScheduleGrid& grid,
552 const FieldPropsManager& fp);
553 void addWell(Well well);
554 void addWell(const std::string& wellName,
555 const std::string& group,
556 int headI,
557 int headJ,
558 Phase preferredPhase,
559 const std::optional<double>& refDepth,
560 double drainageRadius,
561 bool allowCrossFlow,
562 bool automaticShutIn,
563 int pvt_table,
564 WellGasInflowEquation gas_inflow,
565 std::size_t timeStep,
566 ConnectionOrder wellConnectionOrder);
567 bool updateWPAVE(const std::string& wname, std::size_t report_step, const PAvg& pavg);
568
569 void updateGuideRateModel(const GuideRateModel& new_model, std::size_t report_step);
570 GTNode groupTree(const std::string& root_node, std::size_t report_step, std::size_t level, const std::optional<std::string>& parent_name) const;
571 bool checkGroups(const ParseContext& parseContext, ErrorGuard& errors);
572 bool updateWellStatus( const std::string& well, std::size_t reportStep, WellStatus status, std::optional<KeywordLocation> = {});
573 void addWellToGroup( const std::string& group_name, const std::string& well_name , std::size_t timeStep);
574 void iterateScheduleSection(std::size_t load_start,
575 std::size_t load_end,
576 const ParseContext& parseContext,
577 ErrorGuard& errors,
578 const ScheduleGrid& grid,
579 const std::unordered_map<std::string, double> * target_wellpi,
580 const std::string& prefix,
581 const bool keepKeywords,
582 const bool log_to_debug = false);
583 void addACTIONX(const Action::ActionX& action);
584 void addGroupToGroup( const std::string& parent_group, const std::string& child_group);
585 void addGroup(const std::string& groupName , std::size_t timeStep);
586 void addGroup(Group group);
587 void addGroup(const RestartIO::RstGroup& rst_group, std::size_t timeStep);
588 void addWell(const std::string& wellName, const DeckRecord& record,
589 std::size_t timeStep, ConnectionOrder connection_order);
590 void checkIfAllConnectionsIsShut(std::size_t currentStep);
591 void end_report(std::size_t report_step);
594 void handleKeyword(std::size_t currentStep,
595 const ScheduleBlock& block,
596 const DeckKeyword& keyword,
597 const ParseContext& parseContext,
598 ErrorGuard& errors,
599 const ScheduleGrid& grid,
600 const std::vector<std::string>& matching_wells,
601 bool actionx_mode,
602 SimulatorUpdate* sim_update,
603 const std::unordered_map<std::string, double>* target_wellpi,
604 std::unordered_map<std::string, double>& wpimult_global_factor,
605 WelSegsSet* welsegs_wells = nullptr,
606 std::set<std::string>* compsegs_wells = nullptr);
607
608 void internalWELLSTATUSACTIONXFromPYACTION(const std::string& well_name, std::size_t report_step, const std::string& wellStatus);
609 void prefetchPossibleFutureConnections(const ScheduleGrid& grid, const DeckKeyword& keyword);
610 void store_wgnames(const DeckKeyword& keyword);
611 std::vector<std::string> wellNames(const std::string& pattern,
612 const HandlerContext& context,
613 bool allowEmpty = false);
614 std::vector<std::string> wellNames(const std::string& pattern, std::size_t timeStep, const std::vector<std::string>& matching_wells, InputErrorAction error_action, ErrorGuard& errors, const KeywordLocation& location) const;
615 static std::string formatDate(std::time_t t);
616 std::string simulationDays(std::size_t currentStep) const;
617 void applyGlobalWPIMULT( const std::unordered_map<std::string, double>& wpimult_global_factor);
618
619 bool must_write_rst_file(std::size_t report_step) const;
620
621 bool isWList(std::size_t report_step, const std::string& pattern) const;
622
623 SimulatorUpdate applyAction(std::size_t reportStep, const std::string& action_name, const std::vector<std::string>& matching_wells);
624 };
625}
626
627#endif
Definition ActionX.hpp:74
Definition PyAction.hpp:40
Definition State.hpp:40
Definition WGNames.hpp:29
Simple class capturing active cells of a grid.
Definition ActiveGridCells.hpp:36
Definition CompletedCells.hpp:34
Definition DeckKeyword.hpp:36
Definition DeckRecord.hpp:32
Definition Deck.hpp:49
About cell information and dimension: The actual grid information is held in a pointer to an ERT ecl_...
Definition EclipseGrid.hpp:55
Definition EclipseState.hpp:63
Definition ErrorGuard.hpp:29
Definition FieldPropsManager.hpp:42
Definition GTNode.hpp:34
Definition GasLiftOpt.hpp:218
Definition Group.hpp:47
Definition GuideRateModel.hpp:30
Definition HandlerContext.hpp:48
Definition KeywordLocation.hpp:27
Definition PAvg.hpp:30
Definition ParseContext.hpp:84
Definition Runspec.hpp:480
Definition ScheduleBlock.hpp:52
Definition ScheduleDeck.hpp:53
Definition ScheduleGrid.hpp:37
Definition ScheduleState.hpp:94
Definition Schedule.hpp:89
Class for (de-)serializing.
Definition Serializer.hpp:91
bool isSerializing() const
Returns true if we are currently doing a serialization operation.
Definition Serializer.hpp:204
Definition SummaryState.hpp:72
Definition TracerConfig.hpp:33
Definition UDQConfig.hpp:63
Definition UnitSystem.hpp:34
Definition WelSegsSet.hpp:34
Definition WellMatcher.hpp:32
Definition Well.hpp:78
Definition WriteRestartFileEvents.hpp:31
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition group.hpp:33
Definition state.hpp:55
Definition ScheduleStatic.hpp:40
This struct is used to communicate back from the Schedule::applyAction() what needs to be updated in ...
Definition SimulatorUpdate.hpp:33