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