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