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::Properties {
61
62template<class TypeTag, class MyTypeTag>
64 using type = UndefinedProperty;
65};
66
67template <class TypeTag, class MyTypeTag>
69{
70 using type = UndefinedProperty;
71};
72
73template <class TypeTag, class MyTypeTag>
75{
76 using type = UndefinedProperty;
77};
78
79template <class TypeTag, class MyTypeTag>
81{
82 using type = UndefinedProperty;
83};
84
85template <class TypeTag, class MyTypeTag>
87{
88 using type = UndefinedProperty;
89};
90
91template <class TypeTag, class MyTypeTag>
93{
94 using type = UndefinedProperty;
95};
96
97template<class TypeTag>
98struct EnableTerminalOutput<TypeTag, TTag::FlowProblem> {
99 static constexpr bool value = true;
100};
101template<class TypeTag>
103 static constexpr bool value = true;
104};
105
106template <class TypeTag>
108{
109 static constexpr auto* value = "none";
110};
111
112template <class TypeTag>
113struct SaveStep<TypeTag, TTag::FlowProblem>
114{
115 static constexpr auto* value = "";
116};
117
118template <class TypeTag>
119struct SaveFile<TypeTag, TTag::FlowProblem>
120{
121 static constexpr auto* value = "";
122};
123
124template <class TypeTag>
125struct LoadFile<TypeTag, TTag::FlowProblem>
126{
127 static constexpr auto* value = "";
128};
129
130template <class TypeTag>
131struct LoadStep<TypeTag, TTag::FlowProblem>
132{
133 static constexpr int value = -1;
134};
135
136} // namespace Opm::Properties
137
138namespace Opm {
139
141template<class TypeTag>
143{
144public:
145 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
146 using Grid = GetPropType<TypeTag, Properties::Grid>;
147 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
148 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
149 using BlackoilIndices = GetPropType<TypeTag, Properties::Indices>;
150 using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
151 using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
152 using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
153 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
154 using AquiferModel = GetPropType<TypeTag, Properties::AquiferModel>;
155
157 using PolymerModule = BlackOilPolymerModule<TypeTag>;
158 using MICPModule = BlackOilMICPModule<TypeTag>;
159
165
188 : simulator_(simulator)
189 , serializer_(*this,
190 FlowGenericVanguard::comm(),
191 simulator_.vanguard().eclState().getIOConfig(),
192 Parameters::get<TypeTag, Properties::SaveStep>(),
193 Parameters::get<TypeTag, Properties::LoadStep>(),
194 Parameters::get<TypeTag, Properties::SaveFile>(),
195 Parameters::get<TypeTag, Properties::LoadFile>())
196 {
198
199 // Only rank 0 does print to std::cout, and only if specifically requested.
200 this->terminalOutput_ = false;
201 if (this->grid().comm().rank() == 0) {
202 this->terminalOutput_ = Parameters::get<TypeTag, Properties::EnableTerminalOutput>();
203
204 this->startConvergenceOutputThread(Parameters::get<TypeTag, Properties::OutputExtraConvergenceInfo>(),
205 R"(OutputExtraConvergenceInfo (--output-extra-convergence-info))");
206 }
207 }
208
210 {
211 // Safe to call on all ranks, not just the I/O rank.
213 }
214
215 static void registerParameters()
216 {
217 ModelParameters::registerParameters();
218 SolverParameters::registerParameters();
220
221 Parameters::registerParam<TypeTag, Properties::EnableTerminalOutput>
222 ("Print high-level information about the simulation's progress to the terminal");
223 Parameters::registerParam<TypeTag, Properties::EnableAdaptiveTimeStepping>
224 ("Use adaptive time stepping between report steps");
225 Parameters::registerParam<TypeTag, Properties::OutputExtraConvergenceInfo>
226 ("Provide additional convergence output "
227 "files for diagnostic purposes. "
228 "\"none\" gives no extra output and "
229 "overrides all other options, "
230 "\"steps\" generates an INFOSTEP file, "
231 "\"iterations\" generates an INFOITER file. "
232 "Combine options with commas, e.g., "
233 "\"steps,iterations\" for multiple outputs.");
234 Parameters::registerParam<TypeTag, Properties::SaveStep>
235 ("Save serialized state to .OPMRST file. "
236 "Either a specific report step, \"all\" to save "
237 "all report steps or \":x\" to save every x'th step."
238 "Use negative values of \"x\" to keep only the last "
239 "written step, or \"last\" to save every step, keeping "
240 "only the last.");
241 Parameters::registerParam<TypeTag, Properties::LoadStep>
242 ("Load serialized state from .OPMRST file. "
243 "Either a specific report step, or 0 to load last "
244 "stored report step.");
245 Parameters::registerParam<TypeTag, Properties::SaveFile>
246 ("FileName for .OPMRST file used for saving serialized state. "
247 "If empty, CASENAME.OPMRST is used.");
248 Parameters::hideParam<TypeTag, Properties::SaveFile>();
249 Parameters::registerParam<TypeTag, Properties::LoadFile>
250 ("FileName for .OPMRST file used to load serialized state. "
251 "If empty, CASENAME.OPMRST is used.");
252 Parameters::hideParam<TypeTag, Properties::LoadFile>();
253 }
254
262 {
263 init(timer);
264 // Make cache up to date. No need for updating it in elementCtx.
265 simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
266 // Main simulation loop.
267 while (!timer.done()) {
268 bool continue_looping = runStep(timer);
269 if (!continue_looping) break;
270 }
271 return finalize();
272 }
273
274 void init(SimulatorTimer &timer)
275 {
276 simulator_.setEpisodeIndex(-1);
277
278 // Create timers and file for writing timing info.
279 solverTimer_ = std::make_unique<time::StopWatch>();
280 totalTimer_ = std::make_unique<time::StopWatch>();
281 totalTimer_->start();
282
283 // adaptive time stepping
284 bool enableAdaptive = Parameters::get<TypeTag, Properties::EnableAdaptiveTimeStepping>();
285 bool enableTUNING = Parameters::get<TypeTag, Properties::EnableTuning>();
286 if (enableAdaptive) {
287 const UnitSystem& unitSystem = this->simulator_.vanguard().eclState().getUnits();
288 const auto& sched_state = schedule()[timer.currentStepNum()];
289 auto max_next_tstep = sched_state.max_next_tstep(enableTUNING);
290 if (enableTUNING) {
291 adaptiveTimeStepping_ = std::make_unique<TimeStepper>(max_next_tstep,
292 sched_state.tuning(),
293 unitSystem, terminalOutput_);
294 }
295 else {
296 adaptiveTimeStepping_ = std::make_unique<TimeStepper>(unitSystem, max_next_tstep, terminalOutput_);
297 }
298
299 if (isRestart()) {
300 // For restarts the simulator may have gotten some information
301 // about the next timestep size from the OPMEXTRA field
302 adaptiveTimeStepping_->setSuggestedNextStep(simulator_.timeStepSize());
303 }
304 }
305 }
306
307 void updateTUNING(const Tuning& tuning)
308 {
309 modelParam_.tolerance_mb_ = tuning.XXXMBE;
310 if (terminalOutput_) {
311 OpmLog::debug(fmt::format("Setting SimulatorFullyImplicitBlackoil mass balance limit (XXXMBE) to {:.2e}", tuning.XXXMBE));
312 }
313 }
314
316 {
317 if (schedule().exitStatus().has_value()) {
318 if (terminalOutput_) {
319 OpmLog::info("Stopping simulation since EXIT was triggered by an action keyword.");
320 }
321 report_.success.exit_status = schedule().exitStatus().value();
322 return false;
323 }
324
325 if (serializer_.shouldLoad()) {
327 }
328
329 // Report timestep.
330 if (terminalOutput_) {
331 std::ostringstream ss;
332 timer.report(ss);
333 OpmLog::debug(ss.str());
334 }
335
336 if (terminalOutput_) {
338 }
339
340 // write the inital state at the report stage
341 if (timer.initialStep()) {
342 Dune::Timer perfTimer;
343 perfTimer.start();
344
345 simulator_.setEpisodeIndex(-1);
346 simulator_.setEpisodeLength(0.0);
347 simulator_.setTimeStepSize(0.0);
349 simulator_.problem().writeOutput(timer);
350
351 report_.success.output_write_time += perfTimer.stop();
352 }
353
354 // Run a multiple steps of the solver depending on the time step control.
355 solverTimer_->start();
356
357 if (!solver_) {
359 }
360
361 simulator_.startNextEpisode(
362 simulator_.startTime()
363 + schedule().seconds(timer.currentStepNum()),
364 timer.currentStepLength());
365 simulator_.setEpisodeIndex(timer.currentStepNum());
366 if (serializer_.shouldLoad()) {
369 simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
370 }
371 solver_->model().beginReportStep();
372 bool enableTUNING = Parameters::get<TypeTag, Properties::EnableTuning>();
373
374 // If sub stepping is enabled allow the solver to sub cycle
375 // in case the report steps are too large for the solver to converge
376 //
377 // \Note: The report steps are met in any case
378 // \Note: The sub stepping will require a copy of the state variables
380 auto tuningUpdater = [enableTUNING, this, reportStep = timer.currentStepNum()]()
381 {
382 auto& schedule = this->simulator_.vanguard().schedule();
383 auto& events = this->schedule()[reportStep].events();
384
385 if (events.hasEvent(ScheduleEvents::TUNING_CHANGE)) {
386 // Unset the event to not trigger it again on the next sub step
387 schedule.clear_event(ScheduleEvents::TUNING_CHANGE, reportStep);
388 const auto& sched_state = schedule[reportStep];
389 const auto& max_next_tstep = sched_state.max_next_tstep(enableTUNING);
390 const auto& tuning = sched_state.tuning();
391
392 if (enableTUNING) {
393 adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning);
394 // \Note: Assumes TUNING is only used with adaptive time-stepping
395 // \Note: Need to update both solver (model) and simulator since solver is re-created each report step.
396 solver_->model().updateTUNING(tuning);
397 this->updateTUNING(tuning);
398 } else {
399 this->adaptiveTimeStepping_->updateNEXTSTEP(max_next_tstep);
400 }
401 return max_next_tstep >0;
402 }
403 return false;
404 };
405 tuningUpdater();
406 const auto& events = schedule()[timer.currentStepNum()].events();
407 bool event = events.hasEvent(ScheduleEvents::NEW_WELL) ||
408 events.hasEvent(ScheduleEvents::INJECTION_TYPE_CHANGED) ||
409 events.hasEvent(ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER) ||
410 events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE) ||
411 events.hasEvent(ScheduleEvents::INJECTION_UPDATE) ||
412 events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
413 auto stepReport = adaptiveTimeStepping_->step(timer, *solver_, event, nullptr, tuningUpdater);
414 report_ += stepReport;
415 //Pass simulation report to eclwriter for summary output
416 simulator_.problem().setSimulationReport(report_);
417 } else {
418 // solve for complete report step
419 auto stepReport = solver_->step(timer);
420 report_ += stepReport;
421 if (terminalOutput_) {
422 std::ostringstream ss;
423 stepReport.reportStep(ss);
424 OpmLog::info(ss.str());
425 }
426 }
427
428 // write simulation state at the report stage
429 Dune::Timer perfTimer;
430 perfTimer.start();
431 const double nextstep = adaptiveTimeStepping_ ? adaptiveTimeStepping_->suggestedNextStep() : -1.0;
432 simulator_.problem().setNextTimeStepSize(nextstep);
433 simulator_.problem().writeOutput(timer);
434 report_.success.output_write_time += perfTimer.stop();
435
436 solver_->model().endReportStep();
437
438 // take time that was used to solve system for this reportStep
439 solverTimer_->stop();
440
441 // update timing.
442 report_.success.solver_time += solverTimer_->secsSinceStart();
443
444 if (this->grid().comm().rank() == 0) {
445 // Grab the step convergence reports that are new since last we were here.
446 const auto& reps = solver_->model().stepReports();
447 this->writeConvergenceOutput(std::vector<StepReport>{reps.begin() + already_reported_steps_, reps.end()});
448 already_reported_steps_ = reps.size();
449 }
450
451 // Increment timer, remember well state.
452 ++timer;
453
454 if (terminalOutput_) {
455 std::string msg =
456 "Time step took " + std::to_string(solverTimer_->secsSinceStart()) + " seconds; "
457 "total solver time " + std::to_string(report_.success.solver_time) + " seconds.";
458 OpmLog::debug(msg);
459 }
460
461 serializer_.save(timer);
462
463 return true;
464 }
465
467 {
468 // make sure all output is written to disk before run is finished
469 {
470 Dune::Timer finalOutputTimer;
471 finalOutputTimer.start();
472
473 simulator_.problem().finalizeOutput();
474 report_.success.output_write_time += finalOutputTimer.stop();
475 }
476
477 // Stop timer and create timing report
478 totalTimer_->stop();
479 report_.success.total_time = totalTimer_->secsSinceStart();
481
482 return report_;
483 }
484
485 const Grid& grid() const
486 { return simulator_.vanguard().grid(); }
487
488 template<class Serializer>
489 void serializeOp(Serializer& serializer)
490 {
491 serializer(simulator_);
492 serializer(report_);
493 serializer(adaptiveTimeStepping_);
494 }
495
496 const Model& model() const
497 { return solver_->model(); }
498
499protected:
501 void loadState([[maybe_unused]] HDF5Serializer& serializer,
502 [[maybe_unused]] const std::string& groupName) override
503 {
504#if HAVE_HDF5
505 serializer.read(*this, groupName, "simulator_data");
506#endif
507 }
508
510 void saveState([[maybe_unused]] HDF5Serializer& serializer,
511 [[maybe_unused]] const std::string& groupName) const override
512 {
513#if HAVE_HDF5
514 serializer.write(*this, groupName, "simulator_data");
515#endif
516 }
517
519 std::array<std::string,5> getHeader() const override
520 {
521 std::ostringstream str;
522 Parameters::printValues<TypeTag>(str);
523 return {"OPM Flow",
526 simulator_.vanguard().caseName(),
527 str.str()};
528 }
529
531 const std::vector<int>& getCellMapping() const override
532 {
533 return simulator_.vanguard().globalCell();
534 }
535
536
537 std::unique_ptr<Solver> createSolver(WellModel& wellModel)
538 {
539 auto model = std::make_unique<Model>(simulator_,
541 wellModel,
543
544 if (this->modelParam_.write_partitions_) {
545 const auto& iocfg = this->eclState().cfg().io();
546
547 const auto odir = iocfg.getOutputDir()
548 / std::filesystem::path { "partition" }
549 / iocfg.getBaseName();
550
551 if (this->grid().comm().rank() == 0) {
552 create_directories(odir);
553 }
554
555 this->grid().comm().barrier();
556
557 model->writePartitions(odir);
558
559 this->modelParam_.write_partitions_ = false;
560 }
561
562 return std::make_unique<Solver>(solverParam_, std::move(model));
563 }
564
565 const EclipseState& eclState() const
566 { return simulator_.vanguard().eclState(); }
567
568
569 const Schedule& schedule() const
570 { return simulator_.vanguard().schedule(); }
571
572 bool isRestart() const
573 {
574 const auto& initconfig = eclState().getInitConfig();
575 return initconfig.restartRequested();
576 }
577
579 { return simulator_.problem().wellModel(); }
580
581 const WellModel& wellModel_() const
582 { return simulator_.problem().wellModel(); }
583
584 void startConvergenceOutputThread(std::string_view convOutputOptions,
585 std::string_view optionName)
586 {
587 const auto config = ConvergenceOutputConfiguration {
588 convOutputOptions, optionName
589 };
591 return;
592 }
593
595 [compNames = typename Model::ComponentName{}](const int compIdx)
596 { return std::string_view { compNames.name(compIdx) }; }
597 };
598
600 [usys = this->eclState().getUnits()](const double time)
601 { return usys.from_si(UnitSystem::measure::time, time); }
602 };
603
604 this->convergenceOutputQueue_.emplace();
605 this->convergenceOutputObject_.emplace
606 (this->eclState().getIOConfig().getOutputDir(),
607 this->eclState().getIOConfig().getBaseName(),
608 std::move(getPhaseName),
609 std::move(convertTime),
610 config, *this->convergenceOutputQueue_);
611
614 &this->convergenceOutputObject_.value());
615 }
616
617 void writeConvergenceOutput(std::vector<StepReport>&& reports)
618 {
619 if (! this->convergenceOutputThread_.has_value()) {
620 return;
621 }
622
623 auto requests = std::vector<ConvergenceReportQueue::OutputRequest>{};
624 requests.reserve(reports.size());
625
626 for (auto&& report : reports) {
627 requests.push_back({ report.report_step, report.current_step, std::move(report.report) });
628 }
629
630 this->convergenceOutputQueue_->enqueue(std::move(requests));
631 }
632
634 {
635 if (! this->convergenceOutputThread_.has_value()) {
636 return;
637 }
638
639 this->convergenceOutputQueue_->signalLastOutputRequest();
640 this->convergenceOutputThread_->join();
641 }
642
643 // Data.
645
648
649 std::unique_ptr<Solver> solver_;
650
651 // Observed objects.
653 // Misc. data
655
657 std::size_t already_reported_steps_ = 0;
658 std::unique_ptr<time::StopWatch> solverTimer_;
659 std::unique_ptr<time::StopWatch> totalTimer_;
660 std::unique_ptr<TimeStepper> adaptiveTimeStepping_;
661
662 std::optional<ConvergenceReportQueue> convergenceOutputQueue_{};
663 std::optional<ConvergenceOutputThread> convergenceOutputObject_{};
664 std::optional<std::thread> convergenceOutputThread_{};
665
667};
668
669} // namespace Opm
670
671#endif // OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED
Definition: AdaptiveTimeStepping.hpp:223
static void registerParameters()
Definition: AdaptiveTimeStepping.hpp:304
Definition: BlackoilModel.hpp:164
void writePartitions(const std::filesystem::path &odir) const
Definition: BlackoilModel.hpp:1093
BlackoilModelParameters< TypeTag > ModelParameters
Definition: BlackoilModel.hpp:167
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:101
void prepareDeserialize(const int report_step)
Definition: BlackoilWellModel.hpp:252
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:248
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:60
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblem.hpp:112
Class for (de-)serializing using HDF5.
Definition: HDF5Serializer.hpp:37
Definition: NonlinearSolver.hpp:108
a simulator for the blackoil model
Definition: SimulatorFullyImplicitBlackoil.hpp:143
SolverParameters solverParam_
Definition: SimulatorFullyImplicitBlackoil.hpp:647
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: SimulatorFullyImplicitBlackoil.hpp:145
void loadState(HDF5Serializer &serializer, const std::string &groupName) override
Load simulator state from hdf5 serializer.
Definition: SimulatorFullyImplicitBlackoil.hpp:501
void endConvergenceOutputThread()
Definition: SimulatorFullyImplicitBlackoil.hpp:633
PhaseUsage phaseUsage_
Definition: SimulatorFullyImplicitBlackoil.hpp:652
const WellModel & wellModel_() const
Definition: SimulatorFullyImplicitBlackoil.hpp:581
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: SimulatorFullyImplicitBlackoil.hpp:147
std::unique_ptr< TimeStepper > adaptiveTimeStepping_
Definition: SimulatorFullyImplicitBlackoil.hpp:660
GetPropType< TypeTag, Properties::MaterialLawParams > MaterialLawParams
Definition: SimulatorFullyImplicitBlackoil.hpp:153
std::unique_ptr< Solver > createSolver(WellModel &wellModel)
Definition: SimulatorFullyImplicitBlackoil.hpp:537
void startConvergenceOutputThread(std::string_view convOutputOptions, std::string_view optionName)
Definition: SimulatorFullyImplicitBlackoil.hpp:584
bool isRestart() const
Definition: SimulatorFullyImplicitBlackoil.hpp:572
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: SimulatorFullyImplicitBlackoil.hpp:151
typename Model::ModelParameters ModelParameters
Definition: SimulatorFullyImplicitBlackoil.hpp:162
void init(SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:274
void updateTUNING(const Tuning &tuning)
Definition: SimulatorFullyImplicitBlackoil.hpp:307
SimulatorReport run(SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:261
void serializeOp(Serializer &serializer)
Definition: SimulatorFullyImplicitBlackoil.hpp:489
typename Solver::SolverParameters SolverParameters
Definition: SimulatorFullyImplicitBlackoil.hpp:163
SimulatorSerializer serializer_
Definition: SimulatorFullyImplicitBlackoil.hpp:666
const std::vector< int > & getCellMapping() const override
Returns local-to-global cell mapping.
Definition: SimulatorFullyImplicitBlackoil.hpp:531
BlackOilMICPModule< TypeTag > MICPModule
Definition: SimulatorFullyImplicitBlackoil.hpp:158
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: SimulatorFullyImplicitBlackoil.hpp:150
const EclipseState & eclState() const
Definition: SimulatorFullyImplicitBlackoil.hpp:565
std::unique_ptr< time::StopWatch > totalTimer_
Definition: SimulatorFullyImplicitBlackoil.hpp:659
void saveState(HDF5Serializer &serializer, const std::string &groupName) const override
Save simulator state using hdf5 serializer.
Definition: SimulatorFullyImplicitBlackoil.hpp:510
const Grid & grid() const
Definition: SimulatorFullyImplicitBlackoil.hpp:485
void writeConvergenceOutput(std::vector< StepReport > &&reports)
Definition: SimulatorFullyImplicitBlackoil.hpp:617
std::optional< std::thread > convergenceOutputThread_
Definition: SimulatorFullyImplicitBlackoil.hpp:664
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: SimulatorFullyImplicitBlackoil.hpp:152
static void registerParameters()
Definition: SimulatorFullyImplicitBlackoil.hpp:215
GetPropType< TypeTag, Properties::Indices > BlackoilIndices
Definition: SimulatorFullyImplicitBlackoil.hpp:149
const Model & model() const
Definition: SimulatorFullyImplicitBlackoil.hpp:496
SimulatorFullyImplicitBlackoil(Simulator &simulator)
Definition: SimulatorFullyImplicitBlackoil.hpp:187
std::optional< ConvergenceOutputThread > convergenceOutputObject_
Definition: SimulatorFullyImplicitBlackoil.hpp:663
std::unique_ptr< time::StopWatch > solverTimer_
Definition: SimulatorFullyImplicitBlackoil.hpp:658
std::optional< ConvergenceReportQueue > convergenceOutputQueue_
Definition: SimulatorFullyImplicitBlackoil.hpp:662
std::unique_ptr< Solver > solver_
Definition: SimulatorFullyImplicitBlackoil.hpp:649
bool terminalOutput_
Definition: SimulatorFullyImplicitBlackoil.hpp:654
ModelParameters modelParam_
Definition: SimulatorFullyImplicitBlackoil.hpp:646
bool runStep(SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:315
BlackOilPolymerModule< TypeTag > PolymerModule
Definition: SimulatorFullyImplicitBlackoil.hpp:157
GetPropType< TypeTag, Properties::AquiferModel > AquiferModel
Definition: SimulatorFullyImplicitBlackoil.hpp:154
std::size_t already_reported_steps_
Definition: SimulatorFullyImplicitBlackoil.hpp:657
GetPropType< TypeTag, Properties::Grid > Grid
Definition: SimulatorFullyImplicitBlackoil.hpp:146
SimulatorReport finalize()
Definition: SimulatorFullyImplicitBlackoil.hpp:466
SimulatorReport report_
Definition: SimulatorFullyImplicitBlackoil.hpp:656
WellModel & wellModel_()
Definition: SimulatorFullyImplicitBlackoil.hpp:578
const Schedule & schedule() const
Definition: SimulatorFullyImplicitBlackoil.hpp:569
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: SimulatorFullyImplicitBlackoil.hpp:148
~SimulatorFullyImplicitBlackoil()
Definition: SimulatorFullyImplicitBlackoil.hpp:209
std::array< std::string, 5 > getHeader() const override
Returns header data.
Definition: SimulatorFullyImplicitBlackoil.hpp:519
Simulator & simulator_
Definition: SimulatorFullyImplicitBlackoil.hpp:644
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: AluGridVanguard.hpp:57
void outputReportStep(const SimulatorTimer &timer)
Definition: BlackoilPhases.hpp:27
std::string compileTimestamp()
std::string moduleVersion()
PhaseUsage phaseUsageFromDeck(const EclipseState &eclipseState)
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: NonlinearSolver.hpp:114
Definition: BlackoilPhases.hpp:46
Definition: SimulatorFullyImplicitBlackoil.hpp:63
UndefinedProperty type
Definition: SimulatorFullyImplicitBlackoil.hpp:64
Definition: BlackoilWellModel.hpp:84
Definition: SimulatorFullyImplicitBlackoil.hpp:93
UndefinedProperty type
Definition: SimulatorFullyImplicitBlackoil.hpp:94
Definition: SimulatorFullyImplicitBlackoil.hpp:81
UndefinedProperty type
Definition: SimulatorFullyImplicitBlackoil.hpp:82
Definition: SimulatorFullyImplicitBlackoil.hpp:69
UndefinedProperty type
Definition: SimulatorFullyImplicitBlackoil.hpp:70
Definition: SimulatorFullyImplicitBlackoil.hpp:87
UndefinedProperty type
Definition: SimulatorFullyImplicitBlackoil.hpp:88
Definition: SimulatorFullyImplicitBlackoil.hpp:75
UndefinedProperty type
Definition: SimulatorFullyImplicitBlackoil.hpp:76
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