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>
27
28#ifdef RESERVOIR_COUPLING_ENABLED
29#include <opm/input/eclipse/Schedule/ResCoup/ReservoirCouplingInfo.hpp>
30#include <opm/input/eclipse/Schedule/ResCoup/MasterGroup.hpp>
31#include <opm/input/eclipse/Schedule/ResCoup/Slaves.hpp>
34#include <opm/common/Exceptions.hpp>
35#endif
36
37#include <opm/input/eclipse/Units/UnitSystem.hpp>
38
39#include <opm/grid/utility/StopWatch.hpp>
40
54
55#if HAVE_HDF5
57#endif
58
59#include <filesystem>
60#include <memory>
61#include <string>
62#include <string_view>
63#include <utility>
64#include <vector>
65
66#include <fmt/format.h>
67
68namespace Opm::Parameters {
69
70struct EnableAdaptiveTimeStepping { static constexpr bool value = true; };
71struct OutputExtraConvergenceInfo { static constexpr auto* value = "none"; };
72struct SaveStep { static constexpr auto* value = ""; };
73struct SaveFile { static constexpr auto* value = ""; };
74struct LoadFile { static constexpr auto* value = ""; };
75struct LoadStep { static constexpr int value = -1; };
76struct Slave { static constexpr bool value = false; };
77
78} // namespace Opm::Parameters
79
80namespace Opm::detail {
81
83
84}
85
86namespace Opm {
87
89template<class TypeTag>
91{
92protected:
93 struct MPI_Comm_Deleter;
94public:
107
111
113 using ModelParameters = typename Model::ModelParameters;
116
120 : simulator_(simulator)
121 , serializer_(*this,
122 FlowGenericVanguard::comm(),
123 simulator_.vanguard().eclState().getIOConfig(),
124 Parameters::Get<Parameters::SaveStep>(),
125 Parameters::Get<Parameters::LoadStep>(),
126 Parameters::Get<Parameters::SaveFile>(),
127 Parameters::Get<Parameters::LoadFile>())
128 {
129
130 // Only rank 0 does print to std::cout, and only if specifically requested.
131 this->terminalOutput_ = false;
132 if (this->grid().comm().rank() == 0) {
133 this->terminalOutput_ = Parameters::Get<Parameters::EnableTerminalOutput>();
134
136 [compNames = typename Model::ComponentName{}](const int compIdx)
137 { return std::string_view { compNames.name(compIdx) }; }
138 };
139
140 if (!simulator_.vanguard().eclState().getIOConfig().initOnly()) {
142 startThread(this->simulator_.vanguard().eclState(),
143 Parameters::Get<Parameters::OutputExtraConvergenceInfo>(),
144 R"(OutputExtraConvergenceInfo (--output-extra-convergence-info))",
145 getPhaseName);
146 }
147 }
148 }
149
151 {
152 // Safe to call on all ranks, not just the I/O rank.
154 }
155
156 static void registerParameters()
157 {
158 ModelParameters::registerParameters();
159 SolverParameters::registerParameters();
162 }
163
170#ifdef RESERVOIR_COUPLING_ENABLED
171 SimulatorReport run(SimulatorTimer& timer, int argc, char** argv)
172 {
173 init(timer, argc, argv);
174#else
176 {
177 init(timer);
178#endif
179 // Make cache up to date. No need for updating it in elementCtx.
180 // NB! Need to be at the correct step in case of restart
181 simulator_.setEpisodeIndex(timer.currentStepNum());
182 simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
183 // Main simulation loop.
184 while (!timer.done()) {
185 simulator_.problem().writeReports(timer);
186 bool continue_looping = runStep(timer);
187 if (!continue_looping) break;
188 }
189 simulator_.problem().writeReports(timer);
190 return finalize();
191 }
192
193#ifdef RESERVOIR_COUPLING_ENABLED
194 // This method should only be called if slave mode (i.e. Parameters::Get<Parameters::Slave>())
195 // is false. We try to determine if this is a normal flow simulation or a reservoir
196 // coupling master. It is a normal flow simulation if the schedule does not contain
197 // any SLAVES and GRUPMAST keywords.
198 bool checkRunningAsReservoirCouplingMaster()
199 {
200 for (std::size_t report_step = 0; report_step < this->schedule().size(); ++report_step) {
201 auto rescoup = this->schedule()[report_step].rescoup();
202 auto slave_count = rescoup.slaveCount();
203 auto master_group_count = rescoup.masterGroupCount();
204 // - GRUPMAST and SLAVES keywords need to be specified at the same report step
205 // - They can only occur once in the schedule
206 if (slave_count > 0 && master_group_count > 0) {
207 return true;
208 }
209 else if (slave_count > 0 && master_group_count == 0) {
210 throw ReservoirCouplingError(
211 "Inconsistent reservoir coupling master schedule: "
212 "Slave count is greater than 0 but master group count is 0"
213 );
214 }
215 else if (slave_count == 0 && master_group_count > 0) {
216 throw ReservoirCouplingError(
217 "Inconsistent reservoir coupling master schedule: "
218 "Master group count is greater than 0 but slave count is 0"
219 );
220 }
221 }
222 return false;
223 }
224#endif
225
226#ifdef RESERVOIR_COUPLING_ENABLED
227 // NOTE: The argc and argv will be used when launching a slave process
228 void init(const SimulatorTimer& timer, int argc, char** argv)
229 {
230 auto slave_mode = Parameters::Get<Parameters::Slave>();
231 if (slave_mode) {
232 this->reservoirCouplingSlave_ =
233 std::make_unique<ReservoirCouplingSlave<Scalar>>(
235 this->schedule(), timer
236 );
237 this->reservoirCouplingSlave_->sendAndReceiveInitialData();
238 this->simulator_.setReservoirCouplingSlave(this->reservoirCouplingSlave_.get());
239 wellModel_().setReservoirCouplingSlave(this->reservoirCouplingSlave_.get());
240 }
241 else {
242 auto master_mode = checkRunningAsReservoirCouplingMaster();
243 if (master_mode) {
244 this->reservoirCouplingMaster_ =
245 std::make_unique<ReservoirCouplingMaster<Scalar>>(
247 this->schedule(),
248 argc, argv
249 );
250 this->simulator_.setReservoirCouplingMaster(this->reservoirCouplingMaster_.get());
251 wellModel_().setReservoirCouplingMaster(this->reservoirCouplingMaster_.get());
252 }
253 }
254#else
255 void init(const SimulatorTimer& timer)
256 {
257#endif
258 simulator_.setEpisodeIndex(-1);
259
260 // Create timers and file for writing timing info.
261 solverTimer_ = std::make_unique<time::StopWatch>();
262 totalTimer_ = std::make_unique<time::StopWatch>();
263 totalTimer_->start();
264
265 // adaptive time stepping
266 bool enableAdaptive = Parameters::Get<Parameters::EnableAdaptiveTimeStepping>();
267 bool enableTUNING = Parameters::Get<Parameters::EnableTuning>();
268 if (enableAdaptive) {
269 const UnitSystem& unitSystem = this->simulator_.vanguard().eclState().getUnits();
270 const auto& sched_state = schedule()[timer.currentStepNum()];
271 auto max_next_tstep = sched_state.max_next_tstep(enableTUNING);
272 if (enableTUNING) {
273 adaptiveTimeStepping_ = std::make_unique<TimeStepper>(max_next_tstep,
274 sched_state.tuning(),
275 unitSystem, report_, terminalOutput_);
276 }
277 else {
278 adaptiveTimeStepping_ = std::make_unique<TimeStepper>(unitSystem, report_, max_next_tstep, terminalOutput_);
279 }
280 if (isRestart()) {
281 // For restarts the simulator may have gotten some information
282 // about the next timestep size from the OPMEXTRA field
283 adaptiveTimeStepping_->setSuggestedNextStep(simulator_.timeStepSize());
284 }
285 }
286 }
287
288 void updateTUNING(const Tuning& tuning)
289 {
290 modelParam_.tolerance_cnv_ = tuning.TRGCNV;
291 modelParam_.tolerance_cnv_relaxed_ = tuning.XXXCNV;
292 modelParam_.tolerance_mb_ = tuning.TRGMBE;
293 modelParam_.tolerance_mb_relaxed_ = tuning.XXXMBE;
294 modelParam_.newton_max_iter_ = tuning.NEWTMX;
295 modelParam_.newton_min_iter_ = tuning.NEWTMN;
296 if (terminalOutput_) {
297 const auto msg = fmt::format("Tuning values: "
298 "MB: {:.2e}, CNV: {:.2e}, NEWTMN: {}, NEWTMX: {}",
299 tuning.TRGMBE, tuning.TRGCNV, tuning.NEWTMN, tuning.NEWTMX);
300 OpmLog::debug(msg);
301 if (tuning.TRGTTE_has_value) {
302 OpmLog::warning("Tuning item 2-1 (TRGTTE) is not supported.");
303 }
304 if (tuning.TRGLCV_has_value) {
305 OpmLog::warning("Tuning item 2-4 (TRGLCV) is not supported.");
306 }
307 if (tuning.XXXTTE_has_value) {
308 OpmLog::warning("Tuning item 2-5 (XXXTTE) is not supported.");
309 }
310 if (tuning.XXXLCV_has_value) {
311 OpmLog::warning("Tuning item 2-8 (XXXLCV) is not supported.");
312 }
313 if (tuning.XXXWFL_has_value) {
314 OpmLog::warning("Tuning item 2-9 (XXXWFL) is not supported.");
315 }
316 if (tuning.TRGFIP_has_value) {
317 OpmLog::warning("Tuning item 2-10 (TRGFIP) is not supported.");
318 }
319 if (tuning.TRGSFT_has_value) {
320 OpmLog::warning("Tuning item 2-11 (TRGSFT) is not supported.");
321 }
322 if (tuning.THIONX_has_value) {
323 OpmLog::warning("Tuning item 2-12 (THIONX) is not supported.");
324 }
325 if (tuning.TRWGHT_has_value) {
326 OpmLog::warning("Tuning item 2-13 (TRWGHT) is not supported.");
327 }
328 if (tuning.LITMAX_has_value) {
329 OpmLog::warning("Tuning item 3-3 (LITMAX) is not supported.");
330 }
331 if (tuning.LITMIN_has_value) {
332 OpmLog::warning("Tuning item 3-4 (LITMIN) is not supported.");
333 }
334 if (tuning.MXWSIT_has_value) {
335 OpmLog::warning("Tuning item 3-5 (MXWSIT) is not supported.");
336 }
337 if (tuning.MXWPIT_has_value) {
338 OpmLog::warning("Tuning item 3-6 (MXWPIT) is not supported.");
339 }
340 if (tuning.DDPLIM_has_value) {
341 OpmLog::warning("Tuning item 3-7 (DDPLIM) is not supported.");
342 }
343 if (tuning.DDSLIM_has_value) {
344 OpmLog::warning("Tuning item 3-8 (DDSLIM) is not supported.");
345 }
346 if (tuning.TRGDPR_has_value) {
347 OpmLog::warning("Tuning item 3-9 (TRGDPR) is not supported.");
348 }
349 if (tuning.XXXDPR_has_value) {
350 OpmLog::warning("Tuning item 3-10 (XXXDPR) is not supported.");
351 }
352 if (tuning.MNWRFP_has_value) {
353 OpmLog::warning("Tuning item 3-11 (MNWRFP) is not supported.");
354 }
355 }
356 }
357
359 {
360 if (schedule().exitStatus().has_value()) {
361 if (terminalOutput_) {
362 OpmLog::info("Stopping simulation since EXIT was triggered by an action keyword.");
363 }
364 report_.success.exit_status = schedule().exitStatus().value();
365 return false;
366 }
367
368 if (serializer_.shouldLoad()) {
370 }
371
372 // Report timestep.
373 if (terminalOutput_) {
374 std::ostringstream ss;
375 timer.report(ss);
376 OpmLog::debug(ss.str());
378 }
379
380 // write the inital state at the report stage
381 if (timer.initialStep()) {
382 Dune::Timer perfTimer;
383 perfTimer.start();
384
385 simulator_.setEpisodeIndex(-1);
386 simulator_.setEpisodeLength(0.0);
387 simulator_.setTimeStepSize(0.0);
389 simulator_.problem().writeOutput(true);
390
391 report_.success.output_write_time += perfTimer.stop();
392 }
393
394 // Run a multiple steps of the solver depending on the time step control.
395 solverTimer_->start();
396
397 if (!solver_) {
399 }
400
401 simulator_.startNextEpisode(
402 simulator_.startTime()
403 + schedule().seconds(timer.currentStepNum()),
404 timer.currentStepLength());
405 simulator_.setEpisodeIndex(timer.currentStepNum());
406
407 if (serializer_.shouldLoad()) {
410 simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
411 }
412
413 this->solver_->model().beginReportStep();
414
415 const bool enableTUNING = Parameters::Get<Parameters::EnableTuning>();
416
417 // If sub stepping is enabled allow the solver to sub cycle
418 // in case the report steps are too large for the solver to converge
419 //
420 // \Note: The report steps are met in any case
421 // \Note: The sub stepping will require a copy of the state variables
423 auto tuningUpdater = [enableTUNING, this,
424 reportStep = timer.currentStepNum()](const double curr_time,
425 double dt, const int timeStep)
426 {
427 auto& schedule = this->simulator_.vanguard().schedule();
428 auto& events = this->schedule()[reportStep].events();
429
430 bool result = false;
431 if (events.hasEvent(ScheduleEvents::TUNING_CHANGE)) {
432 // Unset the event to not trigger it again on the next sub step
433 schedule.clear_event(ScheduleEvents::TUNING_CHANGE, reportStep);
434 const auto& sched_state = schedule[reportStep];
435 const auto& max_next_tstep = sched_state.max_next_tstep(enableTUNING);
436 const auto& tuning = sched_state.tuning();
437
438 if (enableTUNING) {
439 adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning);
440 // \Note: Assumes TUNING is only used with adaptive time-stepping
441 // \Note: Need to update both solver (model) and simulator since solver is re-created each report step.
442 solver_->model().updateTUNING(tuning);
443 this->updateTUNING(tuning);
444 dt = this->adaptiveTimeStepping_->suggestedNextStep();
445 } else {
446 dt = max_next_tstep;
447 this->adaptiveTimeStepping_->updateNEXTSTEP(max_next_tstep);
448 }
449 result = max_next_tstep > 0;
450 }
451
452 const auto& wcycle = schedule[reportStep].wcycle.get();
453 if (wcycle.empty()) {
454 return result;
455 }
456
457 const auto& wmatcher = schedule.wellMatcher(reportStep);
458 double wcycle_time_step =
459 wcycle.nextTimeStep(curr_time,
460 dt,
461 wmatcher,
462 this->wellModel_().wellOpenTimes(),
463 this->wellModel_().wellCloseTimes(),
464 [timeStep,
465 &wg_events = this->wellModel_().reportStepStartEvents()]
466 (const std::string& name)
467 {
468 if (timeStep != 0) {
469 return false;
470 }
471 return wg_events.hasEvent(name, ScheduleEvents::REQUEST_OPEN_WELL);
472 });
473
474 wcycle_time_step = this->grid().comm().min(wcycle_time_step);
475 if (dt != wcycle_time_step) {
476 this->adaptiveTimeStepping_->updateNEXTSTEP(wcycle_time_step);
477 return true;
478 }
479
480 return result;
481 };
482
483 tuningUpdater(timer.simulationTimeElapsed(),
484 this->adaptiveTimeStepping_->suggestedNextStep(), 0);
485
486#ifdef RESERVOIR_COUPLING_ENABLED
487 if (this->reservoirCouplingMaster_) {
488 this->reservoirCouplingMaster_->maybeSpawnSlaveProcesses(timer.currentStepNum());
489 this->reservoirCouplingMaster_->maybeActivate(timer.currentStepNum());
490 }
491 else if (this->reservoirCouplingSlave_) {
492 this->reservoirCouplingSlave_->maybeActivate(timer.currentStepNum());
493 }
494#endif
495 const auto& events = schedule()[timer.currentStepNum()].events();
496 bool event = events.hasEvent(ScheduleEvents::NEW_WELL) ||
497 events.hasEvent(ScheduleEvents::INJECTION_TYPE_CHANGED) ||
498 events.hasEvent(ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER) ||
499 events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE) ||
500 events.hasEvent(ScheduleEvents::INJECTION_UPDATE) ||
501 events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
502 auto stepReport = adaptiveTimeStepping_->step(timer, *solver_, event, tuningUpdater);
503 report_ += stepReport;
504 //Pass simulation report to eclwriter for summary output
505 simulator_.problem().setSimulationReport(report_);
506 } else {
507 // solve for complete report step
508 auto stepReport = solver_->step(timer, nullptr);
509 report_ += stepReport;
510 if (terminalOutput_) {
511 std::ostringstream ss;
512 stepReport.reportStep(ss);
513 OpmLog::info(ss.str());
514 }
515 }
516
517 // write simulation state at the report stage
518 Dune::Timer perfTimer;
519 perfTimer.start();
520 const double nextstep = adaptiveTimeStepping_ ? adaptiveTimeStepping_->suggestedNextStep() : -1.0;
521 simulator_.problem().setNextTimeStepSize(nextstep);
522 simulator_.problem().writeOutput(true);
523 report_.success.output_write_time += perfTimer.stop();
524
525 solver_->model().endReportStep();
526
527 // take time that was used to solve system for this reportStep
528 solverTimer_->stop();
529
530 // update timing.
531 report_.success.solver_time += solverTimer_->secsSinceStart();
532
533 if (this->grid().comm().rank() == 0) {
534 // Grab the step convergence reports that are new since last we
535 // were here.
536 const auto& reps = this->solver_->model().stepReports();
538 }
539
540 // Increment timer, remember well state.
541 ++timer;
542
543 if (terminalOutput_) {
544 std::string msg =
545 "Time step took " + std::to_string(solverTimer_->secsSinceStart()) + " seconds; "
546 "total solver time " + std::to_string(report_.success.solver_time) + " seconds.";
547 OpmLog::debug(msg);
548 }
549
550 serializer_.save(timer);
551
552 return true;
553 }
554
556 {
557 // make sure all output is written to disk before run is finished
558 {
559 Dune::Timer finalOutputTimer;
560 finalOutputTimer.start();
561
562 simulator_.problem().finalizeOutput();
563 report_.success.output_write_time += finalOutputTimer.stop();
564 }
565
566 // Stop timer and create timing report
567 totalTimer_->stop();
568 report_.success.total_time = totalTimer_->secsSinceStart();
570
571 return report_;
572 }
573
574 const Grid& grid() const
575 { return simulator_.vanguard().grid(); }
576
577 template<class Serializer>
578 void serializeOp(Serializer& serializer)
579 {
580 serializer(simulator_);
581 serializer(report_);
582 serializer(adaptiveTimeStepping_);
583 }
584
585 const Model& model() const
586 { return solver_->model(); }
587
588protected:
590 void loadState([[maybe_unused]] HDF5Serializer& serializer,
591 [[maybe_unused]] const std::string& groupName) override
592 {
593#if HAVE_HDF5
594 serializer.read(*this, groupName, "simulator_data");
595#endif
596 }
597
599 void saveState([[maybe_unused]] HDF5Serializer& serializer,
600 [[maybe_unused]] const std::string& groupName) const override
601 {
602#if HAVE_HDF5
603 serializer.write(*this, groupName, "simulator_data");
604#endif
605 }
606
608 std::array<std::string,5> getHeader() const override
609 {
610 std::ostringstream str;
612 return {"OPM Flow",
615 simulator_.vanguard().caseName(),
616 str.str()};
617 }
618
620 const std::vector<int>& getCellMapping() const override
621 {
622 return simulator_.vanguard().globalCell();
623 }
624
625 std::unique_ptr<Solver> createSolver(WellModel& wellModel)
626 {
627 auto model = std::make_unique<Model>(simulator_,
629 wellModel,
631
632 if (this->modelParam_.write_partitions_) {
633 const auto& iocfg = this->eclState().cfg().io();
634
635 const auto odir = iocfg.getOutputDir()
636 / std::filesystem::path { "partition" }
637 / iocfg.getBaseName();
638
639 if (this->grid().comm().rank() == 0) {
640 create_directories(odir);
641 }
642
643 this->grid().comm().barrier();
644
645 model->writePartitions(odir);
646
647 this->modelParam_.write_partitions_ = false;
648 }
649
650 return std::make_unique<Solver>(solverParam_, std::move(model));
651 }
652
653 const EclipseState& eclState() const
654 { return simulator_.vanguard().eclState(); }
655
656
657 const Schedule& schedule() const
658 { return simulator_.vanguard().schedule(); }
659
660 bool isRestart() const
661 {
662 const auto& initconfig = eclState().getInitConfig();
663 return initconfig.restartRequested();
664 }
665
667 { return simulator_.problem().wellModel(); }
668
669 const WellModel& wellModel_() const
670 { return simulator_.problem().wellModel(); }
671
672 // Data.
674
677
678 std::unique_ptr<Solver> solver_;
679
680 // Observed objects.
681 // Misc. data
683
685 std::unique_ptr<time::StopWatch> solverTimer_;
686 std::unique_ptr<time::StopWatch> totalTimer_;
687 std::unique_ptr<TimeStepper> adaptiveTimeStepping_;
688
690
691#ifdef RESERVOIR_COUPLING_ENABLED
692 bool slaveMode_{false};
693 std::unique_ptr<ReservoirCouplingMaster<Scalar>> reservoirCouplingMaster_{nullptr};
694 std::unique_ptr<ReservoirCouplingSlave<Scalar>> reservoirCouplingSlave_{nullptr};
695#endif
696
698};
699
700} // namespace Opm
701
702#endif // OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED
Definition: AdaptiveTimeStepping.hpp:76
static void registerParameters()
Definition: AdaptiveTimeStepping_impl.hpp:184
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition: blackoilbioeffectsmodules.hh:93
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:64
Class for handling the blackoil well model.
Definition: BlackoilWellModel.hpp:98
void prepareDeserialize(const int report_step)
Definition: BlackoilWellModel.hpp:200
void beginReportStep(const int time_step)
Definition: BlackoilWellModel_impl.hpp:197
std::function< std::string_view(int)> ComponentToPhaseName
Definition: ExtraConvergenceOutputThread.hpp:109
Definition: FlowGenericVanguard.hpp:108
static Parallel::Communication & comm()
Obtain global communicator.
Definition: FlowGenericVanguard.hpp:336
Class for (de-)serializing using HDF5.
Definition: HDF5Serializer.hpp:37
Definition: NonlinearSolver.hpp:96
NonlinearSolverParameters< Scalar > SolverParameters
Definition: NonlinearSolver.hpp:100
Class handling convergence history output for a simulator.
Definition: SimulatorConvergenceOutput.hpp:44
void write(const std::vector< StepReport > &reports)
a simulator for the blackoil model
Definition: SimulatorFullyImplicitBlackoil.hpp:91
SolverParameters solverParam_
Definition: SimulatorFullyImplicitBlackoil.hpp:676
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: SimulatorFullyImplicitBlackoil.hpp:95
void loadState(HDF5Serializer &serializer, const std::string &groupName) override
Load simulator state from hdf5 serializer.
Definition: SimulatorFullyImplicitBlackoil.hpp:590
const WellModel & wellModel_() const
Definition: SimulatorFullyImplicitBlackoil.hpp:669
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: SimulatorFullyImplicitBlackoil.hpp:97
std::unique_ptr< TimeStepper > adaptiveTimeStepping_
Definition: SimulatorFullyImplicitBlackoil.hpp:687
GetPropType< TypeTag, Properties::MaterialLawParams > MaterialLawParams
Definition: SimulatorFullyImplicitBlackoil.hpp:103
std::unique_ptr< Solver > createSolver(WellModel &wellModel)
Definition: SimulatorFullyImplicitBlackoil.hpp:625
void init(const SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:255
bool isRestart() const
Definition: SimulatorFullyImplicitBlackoil.hpp:660
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: SimulatorFullyImplicitBlackoil.hpp:101
typename Model::ModelParameters ModelParameters
Definition: SimulatorFullyImplicitBlackoil.hpp:113
void updateTUNING(const Tuning &tuning)
Definition: SimulatorFullyImplicitBlackoil.hpp:288
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: SimulatorFullyImplicitBlackoil.hpp:106
SimulatorReport run(SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:175
~SimulatorFullyImplicitBlackoil() override
Definition: SimulatorFullyImplicitBlackoil.hpp:150
void serializeOp(Serializer &serializer)
Definition: SimulatorFullyImplicitBlackoil.hpp:578
typename Solver::SolverParameters SolverParameters
Definition: SimulatorFullyImplicitBlackoil.hpp:114
GetPropType< TypeTag, Properties::NonlinearSystem > Model
Definition: SimulatorFullyImplicitBlackoil.hpp:105
SimulatorSerializer serializer_
Definition: SimulatorFullyImplicitBlackoil.hpp:697
const std::vector< int > & getCellMapping() const override
Returns local-to-global cell mapping.
Definition: SimulatorFullyImplicitBlackoil.hpp:620
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: SimulatorFullyImplicitBlackoil.hpp:100
const EclipseState & eclState() const
Definition: SimulatorFullyImplicitBlackoil.hpp:653
std::unique_ptr< time::StopWatch > totalTimer_
Definition: SimulatorFullyImplicitBlackoil.hpp:686
void saveState(HDF5Serializer &serializer, const std::string &groupName) const override
Save simulator state using hdf5 serializer.
Definition: SimulatorFullyImplicitBlackoil.hpp:599
const Grid & grid() const
Definition: SimulatorFullyImplicitBlackoil.hpp:574
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: SimulatorFullyImplicitBlackoil.hpp:102
static void registerParameters()
Definition: SimulatorFullyImplicitBlackoil.hpp:156
GetPropType< TypeTag, Properties::Indices > BlackoilIndices
Definition: SimulatorFullyImplicitBlackoil.hpp:99
const Model & model() const
Definition: SimulatorFullyImplicitBlackoil.hpp:585
SimulatorFullyImplicitBlackoil(Simulator &simulator)
Definition: SimulatorFullyImplicitBlackoil.hpp:119
std::unique_ptr< time::StopWatch > solverTimer_
Definition: SimulatorFullyImplicitBlackoil.hpp:685
SimulatorConvergenceOutput convergence_output_
Definition: SimulatorFullyImplicitBlackoil.hpp:689
std::unique_ptr< Solver > solver_
Definition: SimulatorFullyImplicitBlackoil.hpp:678
bool terminalOutput_
Definition: SimulatorFullyImplicitBlackoil.hpp:682
ModelParameters modelParam_
Definition: SimulatorFullyImplicitBlackoil.hpp:675
bool runStep(SimulatorTimer &timer)
Definition: SimulatorFullyImplicitBlackoil.hpp:358
GetPropType< TypeTag, Properties::AquiferModel > AquiferModel
Definition: SimulatorFullyImplicitBlackoil.hpp:104
GetPropType< TypeTag, Properties::Grid > Grid
Definition: SimulatorFullyImplicitBlackoil.hpp:96
SimulatorReport finalize()
Definition: SimulatorFullyImplicitBlackoil.hpp:555
SimulatorReport report_
Definition: SimulatorFullyImplicitBlackoil.hpp:684
WellModel & wellModel_()
Definition: SimulatorFullyImplicitBlackoil.hpp:666
const Schedule & schedule() const
Definition: SimulatorFullyImplicitBlackoil.hpp:657
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: SimulatorFullyImplicitBlackoil.hpp:98
std::array< std::string, 5 > getHeader() const override
Returns header data.
Definition: SimulatorFullyImplicitBlackoil.hpp:608
Simulator & simulator_
Definition: SimulatorFullyImplicitBlackoil.hpp:673
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:74
void save(SimulatorTimer &timer)
Save data to file if appropriate.
int loadStep() const
Returns step to load.
Definition: SimulatorSerializer.hpp:77
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
double simulationTimeElapsed() const override
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:187
Definition: alignedallocator.hh:32
void registerSimulatorParameters()
void outputReportStep(const SimulatorTimer &timer)
Definition: blackoilbioeffectsmodules.hh:43
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:233
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: SimulatorFullyImplicitBlackoil.hpp:70
static constexpr bool value
Definition: SimulatorFullyImplicitBlackoil.hpp:70
Definition: SimulatorFullyImplicitBlackoil.hpp:74
static constexpr auto * value
Definition: SimulatorFullyImplicitBlackoil.hpp:74
Definition: SimulatorFullyImplicitBlackoil.hpp:75
static constexpr int value
Definition: SimulatorFullyImplicitBlackoil.hpp:75
Definition: SimulatorFullyImplicitBlackoil.hpp:71
static constexpr auto * value
Definition: SimulatorFullyImplicitBlackoil.hpp:71
Definition: SimulatorFullyImplicitBlackoil.hpp:73
static constexpr auto * value
Definition: SimulatorFullyImplicitBlackoil.hpp:73
Definition: SimulatorFullyImplicitBlackoil.hpp:72
static constexpr auto * value
Definition: SimulatorFullyImplicitBlackoil.hpp:72
Definition: SimulatorFullyImplicitBlackoil.hpp:76
static constexpr bool value
Definition: SimulatorFullyImplicitBlackoil.hpp:76
Abstract interface for simulator serialization ops.
Definition: SimulatorSerializer.hpp:36
Definition: SimulatorReport.hpp:122
SimulatorReportSingle success
Definition: SimulatorReport.hpp:123
double solver_time
Definition: SimulatorReport.hpp:38
bool converged
Definition: SimulatorReport.hpp:55
double total_time
Definition: SimulatorReport.hpp:37
int exit_status
Definition: SimulatorReport.hpp:58
double output_write_time
Definition: SimulatorReport.hpp:46