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/ActionResult.hpp>
36#include <opm/input/eclipse/Schedule/Action/SimulatorUpdate.hpp>
37#include <opm/input/eclipse/Schedule/Action/WGNames.hpp>
38#include <opm/input/eclipse/Schedule/CompletedCells.hpp>
39#include <opm/input/eclipse/Schedule/Group/Group.hpp>
40#include <opm/input/eclipse/Schedule/ScheduleDeck.hpp>
41#include <opm/input/eclipse/Schedule/ScheduleState.hpp>
42#include <opm/input/eclipse/Schedule/ScheduleStatic.hpp>
43#include <opm/input/eclipse/Schedule/Well/PAvg.hpp>
44#include <opm/input/eclipse/Schedule/Well/Well.hpp>
45#include <opm/input/eclipse/Schedule/Well/Connection.hpp>
46#include <opm/input/eclipse/Schedule/WriteRestartFileEvents.hpp>
47#include <opm/input/eclipse/Units/UnitSystem.hpp>
48
49namespace Opm {
50 class ActiveGridCells;
51 class Deck;
52 class DeckKeyword;
53 class DeckRecord;
54 enum class ConnectionOrder;
55 class EclipseGrid;
56 class EclipseState;
57 class ErrorGuard;
58 class FieldPropsManager;
59 class GasLiftOpt;
60 class GTNode;
61 class GuideRateConfig;
62 class GuideRateModel;
63 class HandlerContext;
64 enum class InputErrorAction;
65 class NumericalAquifers;
66 class ParseContext;
67 class Python;
68 class Runspec;
69 class RPTConfig;
70 class ScheduleGrid;
71 class SCHEDULESection;
72 class SegmentMatcher;
73 class SummaryState;
74 class TracerConfig;
75 class UDQConfig;
76 class Well;
77 enum class WellGasInflowEquation;
78 class WellMatcher;
79 enum class WellProducerCMode;
80 enum class WellStatus;
81 class WelSegsSet;
82 class WellTestConfig;
83} // namespace Opm
84
85namespace Opm::ReservoirCoupling {
86 class CouplingInfo;
87} // namespace Opm::ReservoirCoupling
88
89namespace Opm::Action {
90 class ActionX;
91 class PyAction;
92 class State;
93} // namespace Opm::Action
94
95namespace Opm::RestartIO {
96 struct RstState;
97} // namespace Opm::RestartIO
98
99namespace Opm {
101 {
102 public:
103 Schedule() = default;
104
105 explicit Schedule(std::shared_ptr<const Python> python_handle);
106
121 Schedule(const Deck& deck,
122 const EclipseGrid& grid,
123 const FieldPropsManager& fp,
124 const NumericalAquifers& numAquifers,
125 const Runspec &runspec,
126 const ParseContext& parseContext,
127 ErrorGuard& errors,
128 std::shared_ptr<const Python> python,
129 const bool lowActionParsingStrictness = false,
130 const bool slave_mode = false,
131 const bool keepKeywords = true,
132 const std::optional<int>& output_interval = {},
133 const RestartIO::RstState* rst = nullptr,
134 const TracerConfig* tracer_config = nullptr);
135
136 template<typename T>
137 Schedule(const Deck& deck,
138 const EclipseGrid& grid,
139 const FieldPropsManager& fp,
140 const NumericalAquifers& numAquifers,
141 const Runspec &runspec,
142 const ParseContext& parseContext,
143 T&& errors,
144 std::shared_ptr<const Python> python,
145 const bool lowActionParsingStrictness = false,
146 const bool slave_mode = false,
147 const bool keepKeywords = true,
148 const std::optional<int>& output_interval = {},
149 const RestartIO::RstState* rst = nullptr,
150 const TracerConfig* tracer_config = nullptr);
151
152 Schedule(const Deck& deck,
153 const EclipseGrid& grid,
154 const FieldPropsManager& fp,
155 const NumericalAquifers& numAquifers,
156 const Runspec &runspec,
157 std::shared_ptr<const Python> python,
158 const bool lowActionParsingStrictness = false,
159 const bool slave_mode = false,
160 const bool keepKeywords = true,
161 const std::optional<int>& output_interval = {},
162 const RestartIO::RstState* rst = nullptr,
163 const TracerConfig* tracer_config = nullptr);
164
165 Schedule(const Deck& deck,
166 const EclipseState& es,
167 const ParseContext& parseContext,
168 ErrorGuard& errors,
169 std::shared_ptr<const Python> python,
170 const bool lowActionParsingStrictness = false,
171 const bool slave_mode = false,
172 const bool keepKeywords = true,
173 const std::optional<int>& output_interval = {},
174 const RestartIO::RstState* rst = nullptr);
175
176 template <typename T>
177 Schedule(const Deck& deck,
178 const EclipseState& es,
179 const ParseContext& parseContext,
180 T&& errors,
181 std::shared_ptr<const Python> python,
182 const bool lowActionParsingStrictness = false,
183 const bool slave_mode = false,
184 const bool keepKeywords = true,
185 const std::optional<int>& output_interval = {},
186 const RestartIO::RstState* rst = nullptr);
187
188 Schedule(const Deck& deck,
189 const EclipseState& es,
190 std::shared_ptr<const Python> python,
191 const bool lowActionParsingStrictness = false,
192 const bool slave_mode = false,
193 const bool keepKeywords = true,
194 const std::optional<int>& output_interval = {},
195 const RestartIO::RstState* rst = nullptr);
196
197 // The constructor *without* the Python arg should really only be used from Python itself
198 Schedule(const Deck& deck,
199 const EclipseState& es,
200 const std::optional<int>& output_interval = {},
201 const RestartIO::RstState* rst = nullptr);
202
203 ~Schedule() = default;
204
205 static Schedule serializationTestObject();
206
207 /*
208 * If the input deck does not specify a start time, Eclipse's 1. Jan
209 * 1983 is defaulted
210 */
211 std::time_t getStartTime() const;
212 std::time_t posixStartTime() const;
213 std::time_t posixEndTime() const;
214 std::time_t simTime(std::size_t timeStep) const;
215 double seconds(std::size_t timeStep) const;
216 double stepLength(std::size_t timeStep) const;
217 std::optional<int> exitStatus() const;
218 const UnitSystem& getUnits() const { return this->m_static.m_unit_system; }
219 const Runspec& runspec() const { return this->m_static.m_runspec; }
220
221 std::size_t numWells() const;
222 std::size_t numWells(std::size_t timestep) const;
223 bool hasWell(const std::string& wellName) const;
224 bool hasWell(const std::string& wellName, std::size_t timeStep) const;
225
226 WellMatcher wellMatcher(std::size_t report_step) const;
227 std::function<std::unique_ptr<SegmentMatcher>()> segmentMatcherFactory(std::size_t report_step) const;
228 std::vector<std::string> wellNames(const std::string& pattern, std::size_t timeStep, const std::vector<std::string>& matching_wells = {}) const;
229 std::vector<std::string> wellNames(const std::string& pattern) const;
230 std::vector<std::string> wellNames(std::size_t timeStep) const;
231 std::vector<std::string> wellNames() const;
240 bool hasGroup(const std::string& groupName, std::size_t timeStep) const;
241
253 std::vector<std::string>
254 groupNames(const std::string& pattern, std::size_t timeStep) const;
255
262 const std::vector<std::string>& groupNames(std::size_t timeStep) const;
263
272 std::vector<std::string> groupNames(const std::string& pattern) const;
273
277 const std::vector<std::string>& groupNames() const;
278
289 std::vector<const Group*> restart_groups(std::size_t timeStep) const;
290
291 std::vector<std::string> changed_wells(std::size_t reportStep) const;
292 const Well& getWell(std::size_t well_index, std::size_t timeStep) const;
293 const Well& getWell(const std::string& wellName, std::size_t timeStep) const;
294 const Well& getWellatEnd(const std::string& well_name) const;
295 // get the list of the constant flux aquifer specified in the whole schedule
296 std::unordered_set<int> getAquiferFluxSchedule() const;
297 std::vector<Well> getWells(std::size_t timeStep) const;
298 std::vector<Well> getWellsatEnd() const;
299
300 // Get wells that have been active any time during simulation
301 std::vector<Well> getActiveWellsAtEnd() const;
302
303 // Get well names of wells that have never been active
304 std::vector<std::string> getInactiveWellNamesAtEnd() const;
305
306 const std::unordered_map<std::string, std::set<int>>&
307 getPossibleFutureConnections() const;
308
309 void shut_well(const std::string& well_name, std::size_t report_step);
310 void shut_well(const std::string& well_name);
311 void stop_well(const std::string& well_name, std::size_t report_step);
312 void stop_well(const std::string& well_name);
313 void open_well(const std::string& well_name, std::size_t report_step);
314 void open_well(const std::string& well_name);
315 void clear_event(ScheduleEvents::Events, std::size_t report_step);
316 void add_event(ScheduleEvents::Events, std::size_t report_step);
317 void applyWellProdIndexScaling(const std::string& well_name, const std::size_t reportStep, const double scalingFactor);
318
320 void clearEvents(const std::size_t report_step);
321
322 WellProducerCMode getGlobalWhistctlMmode(std::size_t timestep) const;
323
324 const UDQConfig& getUDQConfig(std::size_t timeStep) const;
325 void evalAction(const SummaryState& summary_state, std::size_t timeStep);
326
327 GTNode groupTree(std::size_t report_step) const;
328 GTNode groupTree(const std::string& root_node, std::size_t report_step) const;
329 const Group& getGroup(const std::string& groupName, std::size_t timeStep) const;
330
331 std::optional<std::size_t> first_RFT() const;
332 /*
333 Will remove all completions which are connected to cell which is not
334 active. Will scan through all wells and all timesteps.
335 */
336 void filterConnections(const ActiveGridCells& grid);
337 std::size_t size() const;
338
339 bool write_rst_file(std::size_t report_step) const;
340 const std::map< std::string, int >& rst_keywords( size_t timestep ) const;
341
342 // The applyAction() member function is invoked from the simulator
343 // *after* an ACTIONX has triggered. Its return value is a small
344 // structure with 'information' which the simulator should take into
345 // account when updating internal datastructures after the ACTIONX
346 // keywords have been applied.
347 SimulatorUpdate applyAction(std::size_t reportStep,
348 const Action::ActionX& action,
350 const std::unordered_map<std::string, double>& wellpi);
351
352 SimulatorUpdate applyAction(std::size_t reportStep,
353 const Action::ActionX& action,
355 const std::unordered_map<std::string, float>& wellpi);
356 /*
357 The runPyAction() will run the Python script in a PYACTION keyword. In
358 the case of Schedule updates the recommended way of doing that from
359 PYACTION is to invoke a "normal" ACTIONX keyword internally from the
360 Python code. he return value from runPyAction() comes from such a
361 internal ACTIONX.
362 */
363 SimulatorUpdate runPyAction(std::size_t reportStep,
364 const Action::PyAction& pyaction,
365 Action::State& action_state,
366 EclipseState& ecl_state,
367 SummaryState& summary_state);
368
369 SimulatorUpdate runPyAction(std::size_t reportStep,
370 const Action::PyAction& pyaction,
371 Action::State& action_state,
372 EclipseState& ecl_state,
373 SummaryState& summary_state,
374 const std::unordered_map<std::string, double>& target_wellpi);
375
376 SimulatorUpdate runPyAction(std::size_t reportStep,
377 const Action::PyAction& pyaction,
378 Action::State& action_state,
379 EclipseState& ecl_state,
380 SummaryState& summary_state,
381 const std::unordered_map<std::string, float>& target_wellpi);
382
383 SimulatorUpdate modifyCompletions(const std::size_t reportStep,
384 const std::map<std::string, std::vector<Connection>>& extraConns);
385
386 const GasLiftOpt& glo(std::size_t report_step) const;
387
388 bool operator==(const Schedule& data) const;
389 std::shared_ptr<const Python> python() const;
390
391
392 const ScheduleState& back() const;
393 const ScheduleState& operator[](std::size_t index) const;
394 std::vector<ScheduleState>::const_iterator begin() const;
395 std::vector<ScheduleState>::const_iterator end() const;
396 void create_next(const time_point& start_time, const std::optional<time_point>& end_time);
397 void create_next(const ScheduleBlock& block);
398 void create_first(const time_point& start_time, const std::optional<time_point>& end_time);
399
400 void treat_critical_as_non_critical(bool value) { this->m_treat_critical_as_non_critical = value; }
401
402 /*
403 The cmp() function compares two schedule instances in a context aware
404 manner. Floating point numbers are compared with a tolerance. The
405 purpose of this comparison function is to implement regression tests
406 for the schedule instances created by loading a restart file.
407 */
408 static bool cmp(const Schedule& sched1, const Schedule& sched2, std::size_t report_step);
409 void applyKeywords(std::vector<std::unique_ptr<DeckKeyword>>& keywords, std::unordered_map<std::string, double>& target_wellpi, bool action_mode, std::size_t report_step);
410 void applyKeywords(std::vector<std::unique_ptr<DeckKeyword>>& keywords, std::unordered_map<std::string, double>& target_wellpi, bool action_mode);
411
412 template<class Serializer>
413 void serializeOp(Serializer& serializer)
414 {
415 serializer(this->m_static);
416 serializer(this->m_sched_deck);
417 serializer(this->action_wgnames);
418 serializer(this->potential_wellopen_patterns);
419 serializer(this->exit_status);
420 serializer(this->snapshots);
421 serializer(this->restart_output);
422 serializer(this->completed_cells);
423 serializer(this->completed_cells_lgr);
424 serializer(this->completed_cells_lgr_map);
425 serializer(this->m_treat_critical_as_non_critical);
426 serializer(this->current_report_step);
427 serializer(this->m_lowActionParsingStrictness);
428 serializer(this->simUpdateFromPython);
429
430 // If we are deserializing we need to setup the pointer to the
431 // unit system since this is process specific. This is safe
432 // because we set the same value in all well instances.
433 // We do some redundant assignments as these are shared_ptr's
434 // with multiple pointers to any given instance, but it is not
435 // significant so let's keep it simple.
436 if (!serializer.isSerializing()) {
437 for (auto& snapshot : snapshots) {
438 for (auto& well : snapshot.wells) {
439 well.second->updateUnitSystem(&m_static.m_unit_system);
440 }
441 }
442 }
443 }
444
445 template <typename T>
446 std::vector<std::pair<std::size_t, T>> unique() const
447 {
448 std::vector<std::pair<std::size_t, T>> values;
449 for (std::size_t index = 0; index < this->snapshots.size(); index++) {
450 const auto& member = this->snapshots[index].get<T>();
451 const auto& value = member.get();
452 if (values.empty() || !(value == values.back().second))
453 values.push_back( std::make_pair(index, value));
454 }
455 return values;
456 }
457
458 friend std::ostream& operator<<(std::ostream& os, const Schedule& sched);
459 void dump_deck(std::ostream& os) const;
460
461 private:
462 friend class HandlerContext;
463
464 // Please update the member functions
465 // - operator==(const Schedule&) const
466 // - serializationTestObject()
467 // - serializeOp(Serializer&)
468 // when you update/change this list of data members.
469 bool m_treat_critical_as_non_critical = false;
470 ScheduleStatic m_static{};
471 ScheduleDeck m_sched_deck{};
472 Action::WGNames action_wgnames{};
473 std::unordered_set<std::string> potential_wellopen_patterns{}; // Set of well name patterns that potentially can open
474 std::optional<int> exit_status{};
475 std::vector<ScheduleState> snapshots{};
476 WriteRestartFileEvents restart_output{};
477 CompletedCells completed_cells{};
478 std::vector<CompletedCells> completed_cells_lgr{};
479 std::unordered_map<std::string, std::size_t> completed_cells_lgr_map;
480 // Boolean indicating the strictness of parsing process for ActionX and PyAction.
481 // If lowActionParsingStrictness is true, the simulator tries to apply unsupported
482 // keywords, if lowActionParsingStrictness is false, the simulator only applies
483 // supported keywords.
484 bool m_lowActionParsingStrictness = false;
485
486 // This unordered_map contains possible future connections of wells that might get added through an ACTIONX.
487 // For parallel runs, this unordered_map is retrieved by the grid partitioner to ensure these connections
488 // end up on the same partition.
489 std::unordered_map<std::string, std::set<int>> possibleFutureConnections;
490
491 // The current_report_step is set to the current report step when a PYACTION call is executed.
492 // This is needed since the Schedule object does not know the current report step of the simulator and
493 // we only allow PYACTIONS for the current and future report steps.
494 std::size_t current_report_step = 0;
495 // The simUpdateFromPython points to a SimulatorUpdate collecting all updates from one PYACTION call.
496 // The SimulatorUpdate is reset before a new PYACTION call is executed.
497 // It is a shared_ptr, so a Schedule can be constructed using the copy constructor sharing the simUpdateFromPython.
498 // The copy constructor is needed for creating a mocked simulator (msim).
499 std::shared_ptr<SimulatorUpdate> simUpdateFromPython{};
500
501 void init_completed_cells_lgr(const EclipseGrid& ecl_grid);
502 void init_completed_cells_lgr_map(const EclipseGrid& ecl_grid);
503
504 void load_rst(const RestartIO::RstState& rst,
505 const TracerConfig& tracer_config,
506 const ScheduleGrid& grid,
507 const FieldPropsManager& fp);
508 void addWell(Well well);
509 void addWell(const std::string& wellName,
510 const std::string& group,
511 int headI,
512 int headJ,
513 Phase preferredPhase,
514 const std::optional<double>& refDepth,
515 double drainageRadius,
516 bool allowCrossFlow,
517 bool automaticShutIn,
518 int pvt_table,
519 WellGasInflowEquation gas_inflow,
520 std::size_t timeStep,
521 ConnectionOrder wellConnectionOrder);
522 bool updateWPAVE(const std::string& wname, std::size_t report_step, const PAvg& pavg);
523
524 void updateGuideRateModel(const GuideRateModel& new_model, std::size_t report_step);
525 GTNode groupTree(const std::string& root_node, std::size_t report_step, std::size_t level, const std::optional<std::string>& parent_name) const;
526 bool checkGroups(const ParseContext& parseContext, ErrorGuard& errors);
527 bool updateWellStatus( const std::string& well, std::size_t reportStep, WellStatus status, std::optional<KeywordLocation> = {});
528 void addWellToGroup( const std::string& group_name, const std::string& well_name , std::size_t timeStep);
529 void iterateScheduleSection(std::size_t load_start,
530 std::size_t load_end,
531 const ParseContext& parseContext,
532 ErrorGuard& errors,
533 const ScheduleGrid& grid,
534 const std::unordered_map<std::string, double> * target_wellpi,
535 const std::string& prefix,
536 const bool keepKeywords,
537 const bool log_to_debug = false);
538 void addACTIONX(const Action::ActionX& action);
539 void addGroupToGroup( const std::string& parent_group, const std::string& child_group);
540 void addGroup(const std::string& groupName , std::size_t timeStep);
541 void addGroup(Group group);
542 void addGroup(const RestartIO::RstGroup& rst_group, std::size_t timeStep);
543 void addWell(const std::string& wellName, const DeckRecord& record,
544 std::size_t timeStep, ConnectionOrder connection_order);
545 void checkIfAllConnectionsIsShut(std::size_t currentStep);
546 void end_report(std::size_t report_step);
549 void handleKeyword(std::size_t currentStep,
550 const ScheduleBlock& block,
551 const DeckKeyword& keyword,
552 const ParseContext& parseContext,
553 ErrorGuard& errors,
554 const ScheduleGrid& grid,
556 bool action_mode,
557 SimulatorUpdate* sim_update,
558 const std::unordered_map<std::string, double>* target_wellpi,
559 std::unordered_map<std::string, double>& wpimult_global_factor,
560 WelSegsSet* welsegs_wells = nullptr,
561 std::set<std::string>* compsegs_wells = nullptr);
562
563 void internalWELLSTATUSACTIONXFromPYACTION(const std::string& well_name, std::size_t report_step, const std::string& wellStatus);
564 void prefetchPossibleFutureConnections(const ScheduleGrid& grid, const DeckKeyword& keyword,
565 const ParseContext& parseContext, ErrorGuard& errors);
566 void store_wgnames(const DeckKeyword& keyword);
567 std::vector<std::string> wellNames(const std::string& pattern,
568 const HandlerContext& context,
569 bool allowEmpty = false);
570 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;
571 static std::string formatDate(std::time_t t);
572 std::string simulationDays(std::size_t currentStep) const;
573 void applyGlobalWPIMULT( const std::unordered_map<std::string, double>& wpimult_global_factor);
574
575 bool must_write_rst_file(std::size_t report_step) const;
576
577 bool isWList(std::size_t report_step, const std::string& pattern) const;
578
579 SimulatorUpdate applyAction(std::size_t reportStep, const std::string& action_name, const std::vector<std::string>& matching_wells);
580 };
581}
582
583#endif
Definition ActionX.hpp:72
Definition PyAction.hpp:40
Container of matching entities.
Definition ActionResult.hpp:174
Management information about the current run's ACTION system, especially concerning the number of tim...
Definition State.hpp:51
Definition WGNames.hpp:29
Simple class capturing active cells of a grid.
Definition ActiveGridCells.hpp:36
Sparse collection of cells, and their properties, intersected by one or more well connections.
Definition CompletedCells.hpp:36
Definition DeckKeyword.hpp:36
Definition DeckRecord.hpp:32
Definition Deck.hpp:46
About cell information and dimension: The actual grid information is held in a pointer to an ERT ecl_...
Definition EclipseGrid.hpp:62
Definition EclipseState.hpp:63
Definition ErrorGuard.hpp:30
Definition FieldPropsManager.hpp:42
Definition GTNode.hpp:34
Gas lift optimisation parameters for all wells and groups.
Definition GasLiftOpt.hpp:356
Definition Group.hpp:47
Definition GuideRateModel.hpp:30
Definition HandlerContext.hpp:54
Definition KeywordLocation.hpp:27
Definition NumericalAquifers.hpp:38
Definition PAvg.hpp:30
Definition ParseContext.hpp:84
Definition Runspec.hpp:480
Definition ScheduleBlock.hpp:52
Definition ScheduleDeck.hpp:53
Collection of intersected cells and associate properties for all simulation grids,...
Definition ScheduleGrid.hpp:50
Definition ScheduleState.hpp:97
Definition Schedule.hpp:101
bool hasGroup(const std::string &groupName, std::size_t timeStep) const
Query for group existence at particular time.
Definition Schedule.cpp:1188
void clearEvents(const std::size_t report_step)
Clear out all registered events at a given report step.
Definition Schedule.cpp:1003
const std::vector< std::string > & groupNames() const
Retrieve names of all groups in model.
Definition Schedule.cpp:1439
std::vector< const Group * > restart_groups(std::size_t timeStep) const
Retrieve collection of group objects suiteable for restart file output.
Definition Schedule.cpp:1444
Class for (de-)serializing.
Definition Serializer.hpp:94
bool isSerializing() const
Returns true if we are currently doing a serialization operation.
Definition Serializer.hpp:207
Definition SummaryState.hpp:72
Definition TracerConfig.hpp:33
Collection of all user-defined quantities in the current simulation run.
Definition UDQConfig.hpp:69
Definition UnitSystem.hpp:34
Definition WelSegsSet.hpp:34
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:38
Definition state.hpp:56
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