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