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