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#include <opm/simulators/flow/rescoup/ReservoirCouplingEnabled.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>
32#include <opm/simulators/flow/rescoup/ReservoirCouplingMaster.hpp>
33#include <opm/simulators/flow/rescoup/ReservoirCouplingSlave.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
41#include <opm/simulators/aquifers/BlackoilAquiferModel.hpp>
42#include <opm/simulators/flow/BlackoilModel.hpp>
43#include <opm/simulators/flow/BlackoilModelParameters.hpp>
44#include <opm/simulators/flow/ConvergenceOutputConfiguration.hpp>
45#include <opm/simulators/flow/ExtraConvergenceOutputThread.hpp>
46#include <opm/simulators/flow/NonlinearSolver.hpp>
47#include <opm/simulators/flow/SimulatorConvergenceOutput.hpp>
48#include <opm/simulators/flow/SimulatorReportBanners.hpp>
49#include <opm/simulators/flow/SimulatorSerializer.hpp>
50#include <opm/simulators/timestepping/AdaptiveTimeStepping.hpp>
51#include <opm/simulators/timestepping/ConvergenceReport.hpp>
52#include <opm/simulators/utils/moduleVersion.hpp>
53#include <opm/simulators/wells/WellState.hpp>
54
55#if HAVE_HDF5
56#include <opm/simulators/utils/HDF5Serializer.hpp>
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
82void registerSimulatorParameters();
83
84}
85
86namespace Opm {
87
89template<class TypeTag>
91{
92protected:
93 struct MPI_Comm_Deleter;
94public:
99 using BlackoilIndices = GetPropType<TypeTag, Properties::Indices>;
107
108 using TimeStepper = AdaptiveTimeStepping<TypeTag>;
109 using PolymerModule = BlackOilPolymerModule<TypeTag>;
110 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag>;
111
112 using Solver = NonlinearSolver<TypeTag, Model>;
113 using ModelParameters = typename Model::ModelParameters;
114 using SolverParameters = typename Solver::SolverParameters;
115 using WellModel = BlackoilWellModel<TypeTag>;
116
119 explicit SimulatorFullyImplicitBlackoil(Simulator& simulator)
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) {
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()) {
141 this->convergence_output_.
142 startThread(this->simulator_.vanguard().eclState(),
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.
153 convergence_output_.endThread();
154 }
155
156 static void registerParameters()
157 {
158 ModelParameters::registerParameters();
159 SolverParameters::registerParameters();
160 TimeStepper::registerParameters();
161 detail::registerSimulatorParameters();
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
191#ifdef RESERVOIR_COUPLING_ENABLED
192 // Clean up MPI intercommunicators before MPI_Finalize()
193 // Master sends terminate=1 signal; slave receives it and both call MPI_Comm_disconnect()
194 if (this->reservoirCouplingMaster_) {
195 this->reservoirCouplingMaster_->sendTerminateAndDisconnect();
196 }
197 else if (this->reservoirCouplingSlave_ && !this->reservoirCouplingSlave_->terminated()) {
198 // TODO: Implement GECON item 8: stop master process when a slave finishes
199 // Only call if not already terminated via maybeReceiveTerminateSignalFromMaster()
200 // (which happens when master finishes before slave reaches end of its loop)
201 this->reservoirCouplingSlave_->receiveTerminateAndDisconnect();
202 }
203#endif
204
205 return finalize();
206 }
207
208#ifdef RESERVOIR_COUPLING_ENABLED
209 // This method should only be called if slave mode (i.e. Parameters::Get<Parameters::Slave>())
210 // is false. We try to determine if this is a normal flow simulation or a reservoir
211 // coupling master. It is a normal flow simulation if the schedule does not contain
212 // any SLAVES and GRUPMAST keywords.
213 bool checkRunningAsReservoirCouplingMaster()
214 {
215 for (std::size_t report_step = 0; report_step < this->schedule().size(); ++report_step) {
216 auto rescoup = this->schedule()[report_step].rescoup();
217 auto slave_count = rescoup.slaveCount();
218 auto master_group_count = rescoup.masterGroupCount();
219 // Master mode is enabled when SLAVES keyword is present.
220 // - Prediction mode: SLAVES + GRUPMAST (master allocates rates)
221 // - History mode: SLAVES only (master synchronizes time-stepping)
222 if (slave_count > 0) {
223 return true;
224 }
225 else if (master_group_count > 0) {
226 // GRUPMAST without SLAVES is invalid
227 throw ReservoirCouplingError(
228 "Inconsistent reservoir coupling master schedule: "
229 "Master group count is greater than 0 but slave count is 0"
230 );
231 }
232 }
233 return false;
234 }
235#endif
236
237#ifdef RESERVOIR_COUPLING_ENABLED
238 // NOTE: The argc and argv will be used when launching a slave process
239 void init(const SimulatorTimer& timer, int argc, char** argv)
240 {
241 auto slave_mode = Parameters::Get<Parameters::Slave>();
242 if (slave_mode) {
243 this->reservoirCouplingSlave_ =
244 std::make_unique<ReservoirCouplingSlave<Scalar>>(
246 this->schedule(), timer
247 );
248 this->reservoirCouplingSlave_->sendAndReceiveInitialData();
249 this->simulator_.setReservoirCouplingSlave(this->reservoirCouplingSlave_.get());
250 wellModel_().setReservoirCouplingSlave(this->reservoirCouplingSlave_.get());
251 }
252 else {
253 auto master_mode = checkRunningAsReservoirCouplingMaster();
254 if (master_mode) {
255 this->reservoirCouplingMaster_ =
256 std::make_unique<ReservoirCouplingMaster<Scalar>>(
258 this->schedule(),
259 argc, argv
260 );
261 this->simulator_.setReservoirCouplingMaster(this->reservoirCouplingMaster_.get());
262 wellModel_().setReservoirCouplingMaster(this->reservoirCouplingMaster_.get());
263 }
264 }
265#else
266 void init(const SimulatorTimer& timer)
267 {
268#endif
269 simulator_.setEpisodeIndex(-1);
270
271 // Create timers and file for writing timing info.
272 solverTimer_ = std::make_unique<time::StopWatch>();
273 totalTimer_ = std::make_unique<time::StopWatch>();
274 totalTimer_->start();
275
276 // adaptive time stepping
277 bool enableAdaptive = Parameters::Get<Parameters::EnableAdaptiveTimeStepping>();
278 bool enableTUNING = Parameters::Get<Parameters::EnableTuning>();
279 if (enableAdaptive) {
280 const UnitSystem& unitSystem = this->simulator_.vanguard().eclState().getUnits();
281 const auto& sched_state = schedule()[timer.currentStepNum()];
282 auto max_next_tstep = sched_state.max_next_tstep(enableTUNING);
283 if (enableTUNING) {
284 adaptiveTimeStepping_ = std::make_unique<TimeStepper>(max_next_tstep,
285 sched_state.tuning(),
286 unitSystem, report_, terminalOutput_);
287 }
288 else {
289 adaptiveTimeStepping_ = std::make_unique<TimeStepper>(unitSystem, report_, max_next_tstep, terminalOutput_);
290 }
291 if (isRestart()) {
292 // For restarts the simulator may have gotten some information
293 // about the next timestep size from the OPMEXTRA field
294 adaptiveTimeStepping_->setSuggestedNextStep(simulator_.timeStepSize());
295 }
296 }
297 }
298
299 void updateTUNING(const Tuning& tuning)
300 {
301 modelParam_.tolerance_cnv_ = tuning.TRGCNV;
302 modelParam_.tolerance_cnv_relaxed_ = tuning.XXXCNV;
303 modelParam_.tolerance_mb_ = tuning.TRGMBE;
304 modelParam_.tolerance_mb_relaxed_ = tuning.XXXMBE;
305 modelParam_.newton_max_iter_ = tuning.NEWTMX;
306 modelParam_.newton_min_iter_ = tuning.NEWTMN;
307 if (terminalOutput_) {
308 const auto msg = fmt::format(fmt::runtime("Tuning values: "
309 "MB: {:.2e}, CNV: {:.2e}, NEWTMN: {}, NEWTMX: {}"),
310 tuning.TRGMBE, tuning.TRGCNV, tuning.NEWTMN, tuning.NEWTMX);
311 OpmLog::debug(msg);
312 if (tuning.TRGTTE_has_value) {
313 OpmLog::warning("Tuning item 2-1 (TRGTTE) is not supported.");
314 }
315 if (tuning.TRGLCV_has_value) {
316 OpmLog::warning("Tuning item 2-4 (TRGLCV) is not supported.");
317 }
318 if (tuning.XXXTTE_has_value) {
319 OpmLog::warning("Tuning item 2-5 (XXXTTE) is not supported.");
320 }
321 if (tuning.XXXLCV_has_value) {
322 OpmLog::warning("Tuning item 2-8 (XXXLCV) is not supported.");
323 }
324 if (tuning.XXXWFL_has_value) {
325 OpmLog::warning("Tuning item 2-9 (XXXWFL) is not supported.");
326 }
327 if (tuning.TRGFIP_has_value) {
328 OpmLog::warning("Tuning item 2-10 (TRGFIP) is not supported.");
329 }
330 if (tuning.TRGSFT_has_value) {
331 OpmLog::warning("Tuning item 2-11 (TRGSFT) is not supported.");
332 }
333 if (tuning.THIONX_has_value) {
334 OpmLog::warning("Tuning item 2-12 (THIONX) is not supported.");
335 }
336 if (tuning.TRWGHT_has_value) {
337 OpmLog::warning("Tuning item 2-13 (TRWGHT) is not supported.");
338 }
339 if (tuning.LITMAX_has_value) {
340 OpmLog::warning("Tuning item 3-3 (LITMAX) is not supported.");
341 }
342 if (tuning.LITMIN_has_value) {
343 OpmLog::warning("Tuning item 3-4 (LITMIN) is not supported.");
344 }
345 if (tuning.MXWSIT_has_value) {
346 OpmLog::warning("Tuning item 3-5 (MXWSIT) is not supported.");
347 }
348 if (tuning.MXWPIT_has_value) {
349 OpmLog::warning("Tuning item 3-6 (MXWPIT) is not supported.");
350 }
351 if (tuning.DDPLIM_has_value) {
352 OpmLog::warning("Tuning item 3-7 (DDPLIM) is not supported.");
353 }
354 if (tuning.DDSLIM_has_value) {
355 OpmLog::warning("Tuning item 3-8 (DDSLIM) is not supported.");
356 }
357 if (tuning.TRGDPR_has_value) {
358 OpmLog::warning("Tuning item 3-9 (TRGDPR) is not supported.");
359 }
360 if (tuning.XXXDPR_has_value) {
361 OpmLog::warning("Tuning item 3-10 (XXXDPR) is not supported.");
362 }
363 if (tuning.MNWRFP_has_value) {
364 OpmLog::warning("Tuning item 3-11 (MNWRFP) is not supported.");
365 }
366 }
367 }
368
369 void updateTUNINGDP(const TuningDp& tuning_dp)
370 {
371 // NOTE: If TUNINGDP item is _not_ set it should be 0.0
372 modelParam_.tolerance_max_dp_ = tuning_dp.TRGDDP;
373 modelParam_.tolerance_max_ds_ = tuning_dp.TRGDDS;
374 modelParam_.tolerance_max_drs_ = tuning_dp.TRGDDRS;
375 modelParam_.tolerance_max_drv_ = tuning_dp.TRGDDRV;
376
377 // Terminal warnings
378 if (terminalOutput_) {
379 // Warnings unsupported items
380 if (tuning_dp.TRGLCV_has_value) {
381 OpmLog::warning("TUNINGDP item 1 (TRGLCV) is not supported.");
382 }
383 if (tuning_dp.XXXLCV_has_value) {
384 OpmLog::warning("TUNINGDP item 2 (XXXLCV) is not supported.");
385 }
386 }
387 }
388
389 bool runStep(SimulatorTimer& timer)
390 {
391 if (schedule().exitStatus().has_value()) {
392 if (terminalOutput_) {
393 OpmLog::info("Stopping simulation since EXIT was triggered by an action keyword.");
394 }
395 report_.success.exit_status = schedule().exitStatus().value();
396 return false;
397 }
398
399 if (serializer_.shouldLoad()) {
400 serializer_.loadTimerInfo(timer);
401 }
402
403 // Report timestep.
404 if (terminalOutput_) {
405 std::ostringstream ss;
406 timer.report(ss);
407 OpmLog::debug(ss.str());
408 details::outputReportStep(timer);
409 }
410
411 // write the inital state at the report stage
412 if (timer.initialStep()) {
413 Dune::Timer perfTimer;
414 perfTimer.start();
415
416 simulator_.setEpisodeIndex(-1);
417 simulator_.setEpisodeLength(0.0);
418 simulator_.setTimeStepSize(0.0);
419 wellModel_().beginReportStep(timer.currentStepNum());
420 simulator_.problem().writeOutput(true);
421
422 report_.success.output_write_time += perfTimer.stop();
423 }
424
425 // Run a multiple steps of the solver depending on the time step control.
426 solverTimer_->start();
427
428 if (!solver_) {
429 solver_ = createSolver(wellModel_());
430 }
431
432 simulator_.startNextEpisode(
433 simulator_.startTime()
434 + schedule().seconds(timer.currentStepNum()),
435 timer.currentStepLength());
436 simulator_.setEpisodeIndex(timer.currentStepNum());
437
438 if (serializer_.shouldLoad()) {
439 wellModel_().prepareDeserialize(serializer_.loadStep() - 1);
440 serializer_.loadState();
441 simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
442 }
443
444 this->solver_->model().beginReportStep();
445
446 const bool enableTUNING = Parameters::Get<Parameters::EnableTuning>();
447
448 // If sub stepping is enabled allow the solver to sub cycle
449 // in case the report steps are too large for the solver to converge
450 //
451 // \Note: The report steps are met in any case
452 // \Note: The sub stepping will require a copy of the state variables
453 if (adaptiveTimeStepping_) {
454 auto tuningUpdater = [enableTUNING, this,
455 reportStep = timer.currentStepNum()](const double curr_time,
456 double dt, const int timeStep)
457 {
458 auto& schedule = this->simulator_.vanguard().schedule();
459 auto& events = this->schedule()[reportStep].events();
460
461 bool result = false;
462 if (events.hasEvent(ScheduleEvents::TUNING_CHANGE)) {
463 // Unset the event to not trigger it again on the next sub step
464 schedule.clear_event(ScheduleEvents::TUNING_CHANGE, reportStep);
465 const auto& sched_state = schedule[reportStep];
466 const auto& max_next_tstep = sched_state.max_next_tstep(enableTUNING);
467 const auto& tuning = sched_state.tuning();
468
469 if (enableTUNING) {
470 adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning);
471 // \Note: Assumes TUNING is only used with adaptive time-stepping
472 // \Note: Need to update both solver (model) and simulator since solver is re-created each report step.
473 solver_->model().updateTUNING(tuning);
474 this->updateTUNING(tuning);
475 dt = this->adaptiveTimeStepping_->suggestedNextStep();
476 } else {
477 dt = max_next_tstep;
478 this->adaptiveTimeStepping_->updateNEXTSTEP(max_next_tstep);
479 }
480 result = max_next_tstep > 0;
481 }
482
483 if (events.hasEvent(ScheduleEvents::TUNINGDP_CHANGE)) {
484 // Unset the event to not trigger it again on the next sub step
485 schedule.clear_event(ScheduleEvents::TUNINGDP_CHANGE, reportStep);
486
487 // Update TUNINGDP parameters
488 // NOTE: Need to update both solver (model) and simulator since solver is re-created each report
489 // step.
490 const auto& sched_state = schedule[reportStep];
491 const auto& tuning_dp = sched_state.tuning_dp();
492 solver_->model().updateTUNINGDP(tuning_dp);
493 this->updateTUNINGDP(tuning_dp);
494 }
495
496 const auto& wcycle = schedule[reportStep].wcycle.get();
497 if (wcycle.empty()) {
498 return result;
499 }
500
501 const auto& wmatcher = schedule.wellMatcher(reportStep);
502 double wcycle_time_step =
503 wcycle.nextTimeStep(curr_time,
504 dt,
505 wmatcher,
506 this->wellModel_().wellOpenTimes(),
507 this->wellModel_().wellCloseTimes(),
508 [timeStep,
509 &wg_events = this->wellModel_().reportStepStartEvents()]
510 (const std::string& name)
511 {
512 if (timeStep != 0) {
513 return false;
514 }
515 return wg_events.hasEvent(name, ScheduleEvents::REQUEST_OPEN_WELL);
516 });
517
518 wcycle_time_step = this->grid().comm().min(wcycle_time_step);
519 if (dt != wcycle_time_step) {
520 this->adaptiveTimeStepping_->updateNEXTSTEP(wcycle_time_step);
521 return true;
522 }
523
524 return result;
525 };
526
527 tuningUpdater(timer.simulationTimeElapsed(),
528 this->adaptiveTimeStepping_->suggestedNextStep(), 0);
529
530#ifdef RESERVOIR_COUPLING_ENABLED
531 if (this->reservoirCouplingMaster_) {
532 this->reservoirCouplingMaster_->maybeSpawnSlaveProcesses(timer.currentStepNum());
533 this->reservoirCouplingMaster_->maybeActivate(timer.currentStepNum());
534 }
535 else if (this->reservoirCouplingSlave_) {
536 this->reservoirCouplingSlave_->maybeActivate(timer.currentStepNum());
537 }
538#endif
539 const auto& events = schedule()[timer.currentStepNum()].events();
540 bool event = events.hasEvent(ScheduleEvents::NEW_WELL) ||
541 events.hasEvent(ScheduleEvents::INJECTION_TYPE_CHANGED) ||
542 events.hasEvent(ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER) ||
543 events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE) ||
544 events.hasEvent(ScheduleEvents::INJECTION_UPDATE) ||
545 events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
546 auto stepReport = adaptiveTimeStepping_->step(timer, *solver_, event, tuningUpdater);
547 report_ += stepReport;
548 } else {
549 // solve for complete report step
550 auto stepReport = solver_->step(timer, nullptr);
551 report_ += stepReport;
552 // Pass simulation report to eclwriter for summary output
553 simulator_.problem().setSubStepReport(stepReport);
554 simulator_.problem().setSimulationReport(report_);
555 simulator_.problem().endTimeStep();
556 if (terminalOutput_) {
557 std::ostringstream ss;
558 stepReport.reportStep(ss);
559 OpmLog::info(ss.str());
560 }
561 }
562
563 // write simulation state at the report stage
564 Dune::Timer perfTimer;
565 perfTimer.start();
566 const double nextstep = adaptiveTimeStepping_ ? adaptiveTimeStepping_->suggestedNextStep() : -1.0;
567 simulator_.problem().setNextTimeStepSize(nextstep);
568 simulator_.problem().writeOutput(true);
569 report_.success.output_write_time += perfTimer.stop();
570
571 solver_->model().endReportStep();
572
573 // take time that was used to solve system for this reportStep
574 solverTimer_->stop();
575
576 // update timing.
577 report_.success.solver_time += solverTimer_->secsSinceStart();
578
579 if (this->grid().comm().rank() == 0) {
580 // Grab the step convergence reports that are new since last we
581 // were here.
582 const auto& reps = this->solver_->model().stepReports();
583 convergence_output_.write(reps);
584 }
585
586 // Increment timer, remember well state.
587 ++timer;
588
589 if (terminalOutput_) {
590 std::string msg =
591 "Time step took " + std::to_string(solverTimer_->secsSinceStart()) + " seconds; "
592 "total solver time " + std::to_string(report_.success.solver_time) + " seconds.";
593 OpmLog::debug(msg);
594 }
595
596 serializer_.save(timer);
597
598 return true;
599 }
600
601 SimulatorReport finalize()
602 {
603 // make sure all output is written to disk before run is finished
604 {
605 Dune::Timer finalOutputTimer;
606 finalOutputTimer.start();
607
608 simulator_.problem().finalizeOutput();
609 report_.success.output_write_time += finalOutputTimer.stop();
610 }
611
612 // Stop timer and create timing report
613 totalTimer_->stop();
614 report_.success.total_time = totalTimer_->secsSinceStart();
615 report_.success.converged = true;
616
617 return report_;
618 }
619
620 const Grid& grid() const
621 { return simulator_.vanguard().grid(); }
622
623 template<class Serializer>
624 void serializeOp(Serializer& serializer)
625 {
626 serializer(simulator_);
627 serializer(report_);
628 serializer(adaptiveTimeStepping_);
629 }
630
631 const Model& model() const
632 { return solver_->model(); }
633
634protected:
636 void loadState([[maybe_unused]] HDF5Serializer& serializer,
637 [[maybe_unused]] const std::string& groupName) override
638 {
639#if HAVE_HDF5
640 serializer.read(*this, groupName, "simulator_data");
641#endif
642 }
643
645 void saveState([[maybe_unused]] HDF5Serializer& serializer,
646 [[maybe_unused]] const std::string& groupName) const override
647 {
648#if HAVE_HDF5
649 serializer.write(*this, groupName, "simulator_data");
650#endif
651 }
652
654 std::array<std::string,5> getHeader() const override
655 {
656 std::ostringstream str;
658 return {"OPM Flow",
661 simulator_.vanguard().caseName(),
662 str.str()};
663 }
664
666 const std::vector<int>& getCellMapping() const override
667 {
668 return simulator_.vanguard().globalCell();
669 }
670
671 std::unique_ptr<Solver> createSolver(WellModel& wellModel)
672 {
673 auto model = std::make_unique<Model>(simulator_,
674 modelParam_,
675 wellModel,
676 terminalOutput_);
677
678 if (this->modelParam_.write_partitions_) {
679 const auto& iocfg = this->eclState().cfg().io();
680
681 const auto odir = iocfg.getOutputDir()
682 / std::filesystem::path { "partition" }
683 / iocfg.getBaseName();
684
685 if (this->grid().comm().rank() == 0) {
686 create_directories(odir);
687 }
688
689 this->grid().comm().barrier();
690
691 model->writePartitions(odir);
692
693 this->modelParam_.write_partitions_ = false;
694 }
695
696 return std::make_unique<Solver>(solverParam_, std::move(model));
697 }
698
699 const EclipseState& eclState() const
700 { return simulator_.vanguard().eclState(); }
701
702
703 const Schedule& schedule() const
704 { return simulator_.vanguard().schedule(); }
705
706 bool isRestart() const
707 {
708 const auto& initconfig = eclState().getInitConfig();
709 return initconfig.restartRequested();
710 }
711
712 WellModel& wellModel_()
713 { return simulator_.problem().wellModel(); }
714
715 const WellModel& wellModel_() const
716 { return simulator_.problem().wellModel(); }
717
718 // Data.
719 Simulator& simulator_;
720
721 ModelParameters modelParam_;
722 SolverParameters solverParam_;
723
724 std::unique_ptr<Solver> solver_;
725
726 // Observed objects.
727 // Misc. data
728 bool terminalOutput_;
729
730 SimulatorReport report_;
731 std::unique_ptr<time::StopWatch> solverTimer_;
732 std::unique_ptr<time::StopWatch> totalTimer_;
733 std::unique_ptr<TimeStepper> adaptiveTimeStepping_;
734
735 SimulatorConvergenceOutput convergence_output_{};
736
737#ifdef RESERVOIR_COUPLING_ENABLED
738 bool slaveMode_{false};
739 std::unique_ptr<ReservoirCouplingMaster<Scalar>> reservoirCouplingMaster_{nullptr};
740 std::unique_ptr<ReservoirCouplingSlave<Scalar>> reservoirCouplingSlave_{nullptr};
741#endif
742
743 SimulatorSerializer serializer_;
744};
745
746} // namespace Opm
747
748#endif // OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED
Definition AdaptiveTimeStepping.hpp:76
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition blackoilbioeffectsmodules.hh:95
Contains the high level supplements required to extend the black oil model by polymer.
Definition blackoilpolymermodules.hh:65
Class for handling the blackoil well model.
Definition BlackoilWellModel.hpp:98
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:91
void loadState(HDF5Serializer &serializer, const std::string &groupName) override
Load simulator state from hdf5 serializer.
Definition SimulatorFullyImplicitBlackoil.hpp:636
SimulatorReport run(SimulatorTimer &timer)
Run the simulation.
Definition SimulatorFullyImplicitBlackoil.hpp:175
const std::vector< int > & getCellMapping() const override
Returns local-to-global cell mapping.
Definition SimulatorFullyImplicitBlackoil.hpp:666
void saveState(HDF5Serializer &serializer, const std::string &groupName) const override
Save simulator state using hdf5 serializer.
Definition SimulatorFullyImplicitBlackoil.hpp:645
SimulatorFullyImplicitBlackoil(Simulator &simulator)
Initialise from parameters and objects to observe.
Definition SimulatorFullyImplicitBlackoil.hpp:119
std::array< std::string, 5 > getHeader() const override
Returns header data.
Definition SimulatorFullyImplicitBlackoil.hpp:654
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:45
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:742
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187
Definition SimulatorFullyImplicitBlackoil.hpp:70
Definition SimulatorFullyImplicitBlackoil.hpp:74
Definition SimulatorFullyImplicitBlackoil.hpp:75
Definition SimulatorFullyImplicitBlackoil.hpp:71
Definition SimulatorFullyImplicitBlackoil.hpp:73
Definition SimulatorFullyImplicitBlackoil.hpp:72
Definition SimulatorFullyImplicitBlackoil.hpp:76
Abstract interface for simulator serialization ops.
Definition SimulatorSerializer.hpp:36
Definition SimulatorReport.hpp:122