SimulatorFullyImplicitBlackoil.hpp
Go to the documentation of this file.
1/*
2 Copyright 2013, 2015, 2020 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2015 Andreas Lauser
4 Copyright 2017 IRIS
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20*/
21
22#ifndef OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED
23#define OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED
24
25#include <opm/common/ErrorMacros.hpp>
26
27#include <opm/input/eclipse/Units/UnitSystem.hpp>
28
29#include <opm/grid/utility/StopWatch.hpp>
30
43
44#if HAVE_HDF5
46#endif
47
48#include <fmt/format.h>
49
50#include <cstddef>
51#include <filesystem>
52#include <memory>
53#include <optional>
54#include <string>
55#include <string_view>
56#include <thread>
57#include <utility>
58#include <vector>
59
60namespace Opm::Parameters {
61
62struct EnableAdaptiveTimeStepping { static constexpr bool value = true; };
63struct OutputExtraConvergenceInfo { static constexpr auto* value = "none"; };
64struct SaveStep { static constexpr auto* value = ""; };
65struct SaveFile { static constexpr auto* value = ""; };
66struct LoadFile { static constexpr auto* value = ""; };
67struct LoadStep { static constexpr int value = -1; };
68
69} // namespace Opm::Parameters
70
71namespace Opm {
72
74template<class TypeTag>
76{
77public:
88
92
98
121 : simulator_(simulator)
122 , serializer_(*this,
123 FlowGenericVanguard::comm(),
124 simulator_.vanguard().eclState().getIOConfig(),
125 Parameters::Get<Parameters::SaveStep>(),
126 Parameters::Get<Parameters::LoadStep>(),
127 Parameters::Get<Parameters::SaveFile>(),
128 Parameters::Get<Parameters::LoadFile>())
129 {
131
132 // Only rank 0 does print to std::cout, and only if specifically requested.
133 this->terminalOutput_ = false;
134 if (this->grid().comm().rank() == 0) {
135 this->terminalOutput_ = Parameters::Get<Parameters::EnableTerminalOutput>();
136
137 this->startConvergenceOutputThread(Parameters::Get<Parameters::OutputExtraConvergenceInfo>(),
138 R"(OutputExtraConvergenceInfo (--output-extra-convergence-info))");
139 }
140 }
141
143 {
144 // Safe to call on all ranks, not just the I/O rank.
146 }
147
148 static void registerParameters()
149 {
150 ModelParameters::registerParameters();
151 SolverParameters::registerParameters();
153
154 Parameters::Register<Parameters::EnableTerminalOutput>
155 ("Print high-level information about the simulation's progress to the terminal");
156 Parameters::Register<Parameters::EnableAdaptiveTimeStepping>
157 ("Use adaptive time stepping between report steps");
158 Parameters::Register<Parameters::OutputExtraConvergenceInfo>
159 ("Provide additional convergence output "
160 "files for diagnostic purposes. "
161 "\"none\" gives no extra output and "
162 "overrides all other options, "
163 "\"steps\" generates an INFOSTEP file, "
164 "\"iterations\" generates an INFOITER file. "
165 "Combine options with commas, e.g., "
166 "\"steps,iterations\" for multiple outputs.");
167 Parameters::Register<Parameters::SaveStep>
168 ("Save serialized state to .OPMRST file. "
169 "Either a specific report step, \"all\" to save "
170 "all report steps or \":x\" to save every x'th step."
171 "Use negative values of \"x\" to keep only the last "
172 "written step, or \"last\" to save every step, keeping "
173 "only the last.");
174 Parameters::Register<Parameters::LoadStep>
175 ("Load serialized state from .OPMRST file. "
176 "Either a specific report step, or 0 to load last "
177 "stored report step.");
178 Parameters::Register<Parameters::SaveFile>
179 ("FileName for .OPMRST file used for saving serialized state. "
180 "If empty, CASENAME.OPMRST is used.");
181 Parameters::Hide<Parameters::SaveFile>();
182 Parameters::Register<Parameters::LoadFile>
183 ("FileName for .OPMRST file used to load serialized state. "
184 "If empty, CASENAME.OPMRST is used.");
185 Parameters::Hide<Parameters::LoadFile>();
186 }
187
195 {
196 init(timer);
197 // Make cache up to date. No need for updating it in elementCtx.
198 // NB! Need to be at the correct step in case of restart
199 simulator_.setEpisodeIndex(timer.currentStepNum());
200 simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
201 // Main simulation loop.
202 while (!timer.done()) {
203 simulator_.problem().writeReports(timer);
204 bool continue_looping = runStep(timer);
205 if (!continue_looping) break;
206 }
207 simulator_.problem().writeReports(timer);
208 return finalize();
209 }
210
211 void init(SimulatorTimer &timer)
212 {
213 simulator_.setEpisodeIndex(-1);
214
215 // Create timers and file for writing timing info.
216 solverTimer_ = std::make_unique<time::StopWatch>();
217 totalTimer_ = std::make_unique<time::StopWatch>();
218 totalTimer_->start();
219
220 // adaptive time stepping
221 bool enableAdaptive = Parameters::Get<Parameters::EnableAdaptiveTimeStepping>();
222 bool enableTUNING = Parameters::Get<Parameters::EnableTuning>();
223 if (enableAdaptive) {
224 const UnitSystem& unitSystem = this->simulator_.vanguard().eclState().getUnits();
225 const auto& sched_state = schedule()[timer.currentStepNum()];
226 auto max_next_tstep = sched_state.max_next_tstep(enableTUNING);
227 if (enableTUNING) {
228 adaptiveTimeStepping_ = std::make_unique<TimeStepper>(max_next_tstep,
229 sched_state.tuning(),
230 unitSystem, terminalOutput_);
231 }
232 else {
233 adaptiveTimeStepping_ = std::make_unique<TimeStepper>(unitSystem, max_next_tstep, terminalOutput_);
234 }
235
236 if (isRestart()) {
237 // For restarts the simulator may have gotten some information
238 // about the next timestep size from the OPMEXTRA field
239 adaptiveTimeStepping_->setSuggestedNextStep(simulator_.timeStepSize());
240 }
241 }
242 }
243
244 void updateTUNING(const Tuning& tuning)
245 {
246 modelParam_.tolerance_mb_ = tuning.XXXMBE;
247 if (terminalOutput_) {
248 OpmLog::debug(fmt::format("Setting SimulatorFullyImplicitBlackoil mass balance limit (XXXMBE) to {:.2e}", tuning.XXXMBE));
249 }
250 }
251
253 {
254 if (schedule().exitStatus().has_value()) {
255 if (terminalOutput_) {
256 OpmLog::info("Stopping simulation since EXIT was triggered by an action keyword.");
257 }
258 report_.success.exit_status = schedule().exitStatus().value();
259 return false;
260 }
261
262 if (serializer_.shouldLoad()) {
264 }
265
266 // Report timestep.
267 if (terminalOutput_) {
268 std::ostringstream ss;
269 timer.report(ss);
270 OpmLog::debug(ss.str());
271 }
272
273 if (terminalOutput_) {
275 }
276
277 // write the inital state at the report stage
278 if (timer.initialStep()) {
279 Dune::Timer perfTimer;
280 perfTimer.start();
281
282 simulator_.setEpisodeIndex(-1);
283 simulator_.setEpisodeLength(0.0);
284 simulator_.setTimeStepSize(0.0);
286 simulator_.problem().writeOutput();
287
288 report_.success.output_write_time += perfTimer.stop();
289 }
290
291 // Run a multiple steps of the solver depending on the time step control.
292 solverTimer_->start();
293
294 if (!solver_) {
296 }
297
298 simulator_.startNextEpisode(
299 simulator_.startTime()
300 + schedule().seconds(timer.currentStepNum()),
301 timer.currentStepLength());
302 simulator_.setEpisodeIndex(timer.currentStepNum());
303
304 if (serializer_.shouldLoad()) {
307 simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
308 }
309
310 this->solver_->model().beginReportStep();
311
312 const bool enableTUNING = Parameters::Get<Parameters::EnableTuning>();
313
314 // If sub stepping is enabled allow the solver to sub cycle
315 // in case the report steps are too large for the solver to converge
316 //
317 // \Note: The report steps are met in any case
318 // \Note: The sub stepping will require a copy of the state variables
320 auto tuningUpdater = [enableTUNING, this, reportStep = timer.currentStepNum()]()
321 {
322 auto& schedule = this->simulator_.vanguard().schedule();
323 auto& events = this->schedule()[reportStep].events();
324
325 if (events.hasEvent(ScheduleEvents::TUNING_CHANGE)) {
326 // Unset the event to not trigger it again on the next sub step
327 schedule.clear_event(ScheduleEvents::TUNING_CHANGE, reportStep);
328 const auto& sched_state = schedule[reportStep];
329 const auto& max_next_tstep = sched_state.max_next_tstep(enableTUNING);
330 const auto& tuning = sched_state.tuning();
331
332 if (enableTUNING) {
333 adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning);
334 // \Note: Assumes TUNING is only used with adaptive time-stepping
335 // \Note: Need to update both solver (model) and simulator since solver is re-created each report step.
336 solver_->model().updateTUNING(tuning);
337 this->updateTUNING(tuning);
338 } else {
339 this->adaptiveTimeStepping_->updateNEXTSTEP(max_next_tstep);
340 }
341 return max_next_tstep >0;
342 }
343 return false;
344 };
345 tuningUpdater();
346 const auto& events = schedule()[timer.currentStepNum()].events();
347 bool event = events.hasEvent(ScheduleEvents::NEW_WELL) ||
348 events.hasEvent(ScheduleEvents::INJECTION_TYPE_CHANGED) ||
349 events.hasEvent(ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER) ||
350 events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE) ||
351 events.hasEvent(ScheduleEvents::INJECTION_UPDATE) ||
352 events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
353 auto stepReport = adaptiveTimeStepping_->step(timer, *solver_, event, nullptr, tuningUpdater);
354 report_ += stepReport;
355 //Pass simulation report to eclwriter for summary output
356 simulator_.problem().setSimulationReport(report_);
357 } else {
358 // solve for complete report step
359 auto stepReport = solver_->step(timer);
360 report_ += stepReport;
361 if (terminalOutput_) {
362 std::ostringstream ss;
363 stepReport.reportStep(ss);
364 OpmLog::info(ss.str());
365 }
366 }
367
368 // write simulation state at the report stage
369 Dune::Timer perfTimer;
370 perfTimer.start();
371 const double nextstep = adaptiveTimeStepping_ ? adaptiveTimeStepping_->suggestedNextStep() : -1.0;
372 simulator_.problem().setNextTimeStepSize(nextstep);
373 simulator_.problem().writeOutput();
374 report_.success.output_write_time += perfTimer.stop();
375
376 solver_->model().endReportStep();
377
378 // take time that was used to solve system for this reportStep
379 solverTimer_->stop();
380
381 // update timing.
382 report_.success.solver_time += solverTimer_->secsSinceStart();
383
384 if (this->grid().comm().rank() == 0) {
385 // Grab the step convergence reports that are new since last we
386 // were here.
387 const auto& reps = this->solver_->model().stepReports();
388
389 auto reports = std::vector<StepReport> {
390 reps.begin() + this->already_reported_steps_, reps.end()
391 };
392
393 this->writeConvergenceOutput(std::move(reports));
394
395 this->already_reported_steps_ = reps.size();
396 }
397
398 // Increment timer, remember well state.
399 ++timer;
400
401 if (terminalOutput_) {
402 std::string msg =
403 "Time step took " + std::to_string(solverTimer_->secsSinceStart()) + " seconds; "
404 "total solver time " + std::to_string(report_.success.solver_time) + " seconds.";
405 OpmLog::debug(msg);
406 }
407
408 serializer_.save(timer);
409
410 return true;
411 }
412
414 {
415 // make sure all output is written to disk before run is finished
416 {
417 Dune::Timer finalOutputTimer;
418 finalOutputTimer.start();
419
420 simulator_.problem().finalizeOutput();
421 report_.success.output_write_time += finalOutputTimer.stop();
422 }
423
424 // Stop timer and create timing report
425 totalTimer_->stop();
426 report_.success.total_time = totalTimer_->secsSinceStart();
428
429 return report_;
430 }
431
432 const Grid& grid() const
433 { return simulator_.vanguard().grid(); }
434
435 template<class Serializer>
436 void serializeOp(Serializer& serializer)
437 {
438 serializer(simulator_);
439 serializer(report_);
440 serializer(adaptiveTimeStepping_);
441 }
442
443 const Model& model() const
444 { return solver_->model(); }
445
446protected:
448 void loadState([[maybe_unused]] HDF5Serializer& serializer,
449 [[maybe_unused]] const std::string& groupName) override
450 {
451#if HAVE_HDF5
452 serializer.read(*this, groupName, "simulator_data");
453#endif
454 }
455
457 void saveState([[maybe_unused]] HDF5Serializer& serializer,
458 [[maybe_unused]] const std::string& groupName) const override
459 {
460#if HAVE_HDF5
461 serializer.write(*this, groupName, "simulator_data");
462#endif
463 }
464
466 std::array<std::string,5> getHeader() const override
467 {
468 std::ostringstream str;
470 return {"OPM Flow",
473 simulator_.vanguard().caseName(),
474 str.str()};
475 }
476
478 const std::vector<int>& getCellMapping() const override
479 {
480 return simulator_.vanguard().globalCell();
481 }
482
483
484 std::unique_ptr<Solver> createSolver(WellModel& wellModel)
485 {
486 auto model = std::make_unique<Model>(simulator_,
488 wellModel,
490
491 if (this->modelParam_.write_partitions_) {
492 const auto& iocfg = this->eclState().cfg().io();
493
494 const auto odir = iocfg.getOutputDir()
495 / std::filesystem::path { "partition" }
496 / iocfg.getBaseName();
497
498 if (this->grid().comm().rank() == 0) {
499 create_directories(odir);
500 }
501
502 this->grid().comm().barrier();
503
504 model->writePartitions(odir);
505
506 this->modelParam_.write_partitions_ = false;
507 }
508
509 return std::make_unique<Solver>(solverParam_, std::move(model));
510 }
511
512 const EclipseState& eclState() const
513 { return simulator_.vanguard().eclState(); }
514
515
516 const Schedule& schedule() const
517 { return simulator_.vanguard().schedule(); }
518
519 bool isRestart() const
520 {
521 const auto& initconfig = eclState().getInitConfig();
522 return initconfig.restartRequested();
523 }
524
526 { return simulator_.problem().wellModel(); }
527
528 const WellModel& wellModel_() const
529 { return simulator_.problem().wellModel(); }
530
531 void startConvergenceOutputThread(std::string_view convOutputOptions,
532 std::string_view optionName)
533 {
534 const auto config = ConvergenceOutputConfiguration {
535 convOutputOptions, optionName
536 };
538 return;
539 }
540
542 [compNames = typename Model::ComponentName{}](const int compIdx)
543 { return std::string_view { compNames.name(compIdx) }; }
544 };
545
547 [usys = this->eclState().getUnits()](const double time)
548 { return usys.from_si(UnitSystem::measure::time, time); }
549 };
550
551 this->convergenceOutputQueue_.emplace();
552 this->convergenceOutputObject_.emplace
553 (this->eclState().getIOConfig().getOutputDir(),
554 this->eclState().getIOConfig().getBaseName(),
555 std::move(getPhaseName),
556 std::move(convertTime),
557 config, *this->convergenceOutputQueue_);
558
561 &this->convergenceOutputObject_.value());
562 }
563
564 void writeConvergenceOutput(std::vector<StepReport>&& reports)
565 {
566 if (! this->convergenceOutputThread_.has_value()) {
567 return;
568 }
569
570 auto requests = std::vector<ConvergenceReportQueue::OutputRequest>{};
571 requests.reserve(reports.size());
572
573 for (auto&& report : reports) {
574 requests.push_back({ report.report_step, report.current_step, std::move(report.report) });
575 }
576
577 this->convergenceOutputQueue_->enqueue(std::move(requests));
578 }
579
581 {
582 if (! this->convergenceOutputThread_.has_value()) {
583 return;
584 }
585
586 this->convergenceOutputQueue_->signalLastOutputRequest();
587 this->convergenceOutputThread_->join();
588 }
589
590 // Data.
592
595
596 std::unique_ptr<Solver> solver_;
597
598 // Observed objects.
600 // Misc. data
602
604 std::size_t already_reported_steps_ = 0;
605 std::unique_ptr<time::StopWatch> solverTimer_;
606 std::unique_ptr<time::StopWatch> totalTimer_;
607 std::unique_ptr<TimeStepper> adaptiveTimeStepping_;
608
609 std::optional<ConvergenceReportQueue> convergenceOutputQueue_{};
610 std::optional<ConvergenceOutputThread> convergenceOutputObject_{};
611 std::optional<std::thread> convergenceOutputThread_{};
612
614};
615
616} // namespace Opm
617
618#endif // OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED
Definition: AdaptiveTimeStepping.hpp:90
static void registerParameters()
Definition: AdaptiveTimeStepping.hpp:172
Contains the high level supplements required to extend the black oil model by MICP.
Definition: blackoilmicpmodules.hh:49
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:54
Definition: BlackoilModel.hpp:163
BlackoilModelParameters< Scalar > ModelParameters
Definition: BlackoilModel.hpp:178
void writePartitions(const std::filesystem::path &odir) const
Definition: BlackoilModel.hpp:1183
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:98
void prepareDeserialize(const int report_step)
Definition: BlackoilWellModel.hpp:250
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:274
Definition: ComponentName.hpp:34
Definition: ConvergenceOutputConfiguration.hpp:46
std::function< std::string_view(int)> ComponentToPhaseName
Definition: ExtraConvergenceOutputThread.hpp:109
std::function< double(double)> ConvertToTimeUnits
Definition: ExtraConvergenceOutputThread.hpp:115
Definition: FlowGenericVanguard.hpp:99
Class for (de-)serializing using HDF5.
Definition: HDF5Serializer.hpp:37
Definition: NonlinearSolver.hpp:84
a simulator for the blackoil model
Definition: SimulatorFullyImplicitBlackoil.hpp:76
SolverParameters solverParam_
Definition: SimulatorFullyImplicitBlackoil.hpp:594
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: SimulatorFullyImplicitBlackoil.hpp:78
void loadState(HDF5Serializer &serializer, const std::string &groupName) override
Load simulator state from hdf5 serializer.
Definition: SimulatorFullyImplicitBlackoil.hpp:448
void endConvergenceOutputThread()
Definition: SimulatorFullyImplicitBlackoil.hpp:580
PhaseUsage phaseUsage_
Definition: SimulatorFullyImplicitBlackoil.hpp:599
const WellModel & wellModel_() const
Definition: SimulatorFullyImplicitBlackoil.hpp:528
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: SimulatorFullyImplicitBlackoil.hpp:80
std::unique_ptr< TimeStepper > adaptiveTimeStepping_
Definition: SimulatorFullyImplicitBlackoil.hpp:607
GetPropType< TypeTag, Properties::MaterialLawParams > MaterialLawParams
Definition: SimulatorFullyImplicitBlackoil.hpp:86
std::unique_ptr< Solver > createSolver(WellModel &wellModel)
Definition: SimulatorFullyImplicitBlackoil.hpp:484
void startConvergenceOutputThread(std::string_view convOutputOptions, std::string_view optionName)
Definition: SimulatorFullyImplicitBlackoil.hpp:531
bool isRestart() const
Definition: SimulatorFullyImplicitBlackoil.hpp:519
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: SimulatorFullyImplicitBlackoil.hpp:84
typename Model::ModelParameters ModelParameters
Definition: SimulatorFullyImplicitBlackoil.hpp:95
void init(SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:211
void updateTUNING(const Tuning &tuning)
Definition: SimulatorFullyImplicitBlackoil.hpp:244
SimulatorReport run(SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:194
void serializeOp(Serializer &serializer)
Definition: SimulatorFullyImplicitBlackoil.hpp:436
typename Solver::SolverParameters SolverParameters
Definition: SimulatorFullyImplicitBlackoil.hpp:96
SimulatorSerializer serializer_
Definition: SimulatorFullyImplicitBlackoil.hpp:613
const std::vector< int > & getCellMapping() const override
Returns local-to-global cell mapping.
Definition: SimulatorFullyImplicitBlackoil.hpp:478
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: SimulatorFullyImplicitBlackoil.hpp:83
const EclipseState & eclState() const
Definition: SimulatorFullyImplicitBlackoil.hpp:512
std::unique_ptr< time::StopWatch > totalTimer_
Definition: SimulatorFullyImplicitBlackoil.hpp:606
void saveState(HDF5Serializer &serializer, const std::string &groupName) const override
Save simulator state using hdf5 serializer.
Definition: SimulatorFullyImplicitBlackoil.hpp:457
const Grid & grid() const
Definition: SimulatorFullyImplicitBlackoil.hpp:432
void writeConvergenceOutput(std::vector< StepReport > &&reports)
Definition: SimulatorFullyImplicitBlackoil.hpp:564
std::optional< std::thread > convergenceOutputThread_
Definition: SimulatorFullyImplicitBlackoil.hpp:611
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: SimulatorFullyImplicitBlackoil.hpp:85
static void registerParameters()
Definition: SimulatorFullyImplicitBlackoil.hpp:148
GetPropType< TypeTag, Properties::Indices > BlackoilIndices
Definition: SimulatorFullyImplicitBlackoil.hpp:82
const Model & model() const
Definition: SimulatorFullyImplicitBlackoil.hpp:443
SimulatorFullyImplicitBlackoil(Simulator &simulator)
Definition: SimulatorFullyImplicitBlackoil.hpp:120
std::optional< ConvergenceOutputThread > convergenceOutputObject_
Definition: SimulatorFullyImplicitBlackoil.hpp:610
std::unique_ptr< time::StopWatch > solverTimer_
Definition: SimulatorFullyImplicitBlackoil.hpp:605
std::optional< ConvergenceReportQueue > convergenceOutputQueue_
Definition: SimulatorFullyImplicitBlackoil.hpp:609
std::unique_ptr< Solver > solver_
Definition: SimulatorFullyImplicitBlackoil.hpp:596
bool terminalOutput_
Definition: SimulatorFullyImplicitBlackoil.hpp:601
ModelParameters modelParam_
Definition: SimulatorFullyImplicitBlackoil.hpp:593
bool runStep(SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:252
GetPropType< TypeTag, Properties::AquiferModel > AquiferModel
Definition: SimulatorFullyImplicitBlackoil.hpp:87
std::size_t already_reported_steps_
Definition: SimulatorFullyImplicitBlackoil.hpp:604
GetPropType< TypeTag, Properties::Grid > Grid
Definition: SimulatorFullyImplicitBlackoil.hpp:79
SimulatorReport finalize()
Definition: SimulatorFullyImplicitBlackoil.hpp:413
SimulatorReport report_
Definition: SimulatorFullyImplicitBlackoil.hpp:603
WellModel & wellModel_()
Definition: SimulatorFullyImplicitBlackoil.hpp:525
const Schedule & schedule() const
Definition: SimulatorFullyImplicitBlackoil.hpp:516
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: SimulatorFullyImplicitBlackoil.hpp:81
~SimulatorFullyImplicitBlackoil()
Definition: SimulatorFullyImplicitBlackoil.hpp:142
std::array< std::string, 5 > getHeader() const override
Returns header data.
Definition: SimulatorFullyImplicitBlackoil.hpp:466
Simulator & simulator_
Definition: SimulatorFullyImplicitBlackoil.hpp:591
Class handling simulator serialization.
Definition: SimulatorSerializer.hpp:55
void loadState()
Load state from file.
void loadTimerInfo(SimulatorTimer &timer)
Loads time step info from file.
bool shouldLoad() const
Returns whether or not a state should be loaded.
Definition: SimulatorSerializer.hpp:71
void save(SimulatorTimer &timer)
Save data to file if appropriate.
int loadStep() const
Returns step to load.
Definition: SimulatorSerializer.hpp:74
Definition: SimulatorTimer.hpp:39
double currentStepLength() const override
bool initialStep() const override
Whether the current step is the first step.
void report(std::ostream &os) const
int currentStepNum() const override
bool done() const override
Return true if op++() has been called numSteps() times.
Definition: blackoilnewtonmethodparams.hpp:31
void printValues(std::ostream &os)
Print values of the run-time parameters.
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:185
void outputReportStep(const SimulatorTimer &timer)
Definition: blackoilboundaryratevector.hh:37
std::string compileTimestamp()
std::string moduleVersion()
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
PhaseUsage phaseUsageFromDeck(const EclipseState &eclipseState)
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: NonlinearSolver.hpp:90
Definition: SimulatorFullyImplicitBlackoil.hpp:62
static constexpr bool value
Definition: SimulatorFullyImplicitBlackoil.hpp:62
Definition: SimulatorFullyImplicitBlackoil.hpp:66
static constexpr auto * value
Definition: SimulatorFullyImplicitBlackoil.hpp:66
Definition: SimulatorFullyImplicitBlackoil.hpp:67
static constexpr int value
Definition: SimulatorFullyImplicitBlackoil.hpp:67
Definition: SimulatorFullyImplicitBlackoil.hpp:63
static constexpr auto * value
Definition: SimulatorFullyImplicitBlackoil.hpp:63
Definition: SimulatorFullyImplicitBlackoil.hpp:65
static constexpr auto * value
Definition: SimulatorFullyImplicitBlackoil.hpp:65
Definition: SimulatorFullyImplicitBlackoil.hpp:64
static constexpr auto * value
Definition: SimulatorFullyImplicitBlackoil.hpp:64
Definition: BlackoilPhases.hpp:46
Abstract interface for simulator serialization ops.
Definition: SimulatorSerializer.hpp:36
Definition: SimulatorReport.hpp:100
SimulatorReportSingle success
Definition: SimulatorReport.hpp:101
double solver_time
Definition: SimulatorReport.hpp:38
bool converged
Definition: SimulatorReport.hpp:54
double total_time
Definition: SimulatorReport.hpp:37
int exit_status
Definition: SimulatorReport.hpp:56
double output_write_time
Definition: SimulatorReport.hpp:45