/* Copyright 2013, 2015, 2020 SINTEF Digital, Mathematics and Cybernetics. Copyright 2015 Andreas Lauser Copyright 2017 IRIS This file is part of the Open Porous Media project (OPM). OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OPM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OPM. If not, see . */ #ifndef OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED #define OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #if HAVE_HDF5 #include #endif #include #include #include #include #include #include #include namespace Opm::Parameters { struct EnableAdaptiveTimeStepping { static constexpr bool value = true; }; struct OutputExtraConvergenceInfo { static constexpr auto* value = "none"; }; struct SaveStep { static constexpr auto* value = ""; }; struct SaveFile { static constexpr auto* value = ""; }; struct LoadFile { static constexpr auto* value = ""; }; struct LoadStep { static constexpr int value = -1; }; struct Slave { static constexpr bool value = false; }; } // namespace Opm::Parameters namespace Opm::detail { void registerSimulatorParameters(); } namespace Opm { /// a simulator for the blackoil model template class SimulatorFullyImplicitBlackoil : private SerializableSim { protected: struct MPI_Comm_Deleter; public: using Simulator = GetPropType; using Grid = GetPropType; using FluidSystem = GetPropType; using ElementContext = GetPropType; using BlackoilIndices = GetPropType; using PrimaryVariables = GetPropType; using MaterialLaw = GetPropType; using SolutionVector = GetPropType; using MaterialLawParams = GetPropType; using AquiferModel = GetPropType; using TimeStepper = AdaptiveTimeStepping; using PolymerModule = BlackOilPolymerModule; using MICPModule = BlackOilMICPModule; using Model = BlackoilModel; using Solver = NonlinearSolver; using ModelParameters = typename Model::ModelParameters; using SolverParameters = typename Solver::SolverParameters; using WellModel = BlackoilWellModel; /// Initialise from parameters and objects to observe. /// \param[in] param parameters, this class accepts the following: /// parameter (default) effect /// ----------------------------------------------------------- /// output (true) write output to files? /// output_dir ("output") output directoty /// output_interval (1) output every nth step /// nl_pressure_residual_tolerance (0.0) pressure solver residual tolerance (in Pascal) /// nl_pressure_change_tolerance (1.0) pressure solver change tolerance (in Pascal) /// nl_pressure_maxiter (10) max nonlinear iterations in pressure /// nl_maxiter (30) max nonlinear iterations in transport /// nl_tolerance (1e-9) transport solver absolute residual tolerance /// num_transport_substeps (1) number of transport steps per pressure step /// use_segregation_split (false) solve for gravity segregation (if false, /// segregation is ignored). /// /// \param[in] props fluid and rock properties /// \param[in] linsolver linear solver /// \param[in] eclipse_state the object which represents an internalized ECL deck /// \param[in] output_writer /// \param[in] threshold_pressures_by_face if nonempty, threshold pressures that inhibit flow explicit SimulatorFullyImplicitBlackoil(Simulator& simulator) : simulator_(simulator) , convergence_output_(simulator_.vanguard().eclState()) , serializer_(*this, FlowGenericVanguard::comm(), simulator_.vanguard().eclState().getIOConfig(), Parameters::Get(), Parameters::Get(), Parameters::Get(), Parameters::Get()) { phaseUsage_ = phaseUsageFromDeck(eclState()); // Only rank 0 does print to std::cout, and only if specifically requested. this->terminalOutput_ = false; if (this->grid().comm().rank() == 0) { this->terminalOutput_ = Parameters::Get(); auto getPhaseName = ConvergenceOutputThread::ComponentToPhaseName { [compNames = typename Model::ComponentName{}](const int compIdx) { return std::string_view { compNames.name(compIdx) }; } }; convergence_output_. startThread(Parameters::Get(), R"(OutputExtraConvergenceInfo (--output-extra-convergence-info))", getPhaseName); } } ~SimulatorFullyImplicitBlackoil() { // Safe to call on all ranks, not just the I/O rank. convergence_output_.endThread(); } static void registerParameters() { ModelParameters::registerParameters(); SolverParameters::registerParameters(); TimeStepper::registerParameters(); detail::registerSimulatorParameters(); } /// Run the simulation. /// This will run succesive timesteps until timer.done() is true. It will /// modify the reservoir and well states. /// \param[in,out] timer governs the requested reporting timesteps /// \param[in,out] state state of reservoir: pressure, fluxes /// \return simulation report, with timing data SimulatorReport run(SimulatorTimer& timer, int argc, char** argv) { init(timer, argc, argv); // Make cache up to date. No need for updating it in elementCtx. // NB! Need to be at the correct step in case of restart simulator_.setEpisodeIndex(timer.currentStepNum()); simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0); // Main simulation loop. while (!timer.done()) { simulator_.problem().writeReports(timer); bool continue_looping = runStep(timer); if (!continue_looping) break; } simulator_.problem().writeReports(timer); return finalize(); } // NOTE: The argc and argv will be used when launching a slave process void init(SimulatorTimer &timer, int argc, char** argv) { auto slave_mode = Parameters::get(); if (slave_mode) { this->reservoirCouplingSlave_ = std::make_unique( FlowGenericVanguard::comm(), this->schedule() ); this->reservoirCouplingSlave_->sendSimulationStartDateToMasterProcess(); } else { // For now, we require that SLAVES and GRUPMAST are defined at the first // schedule step, so it is enough to check the first step. See the // keyword handlers in opm-common for more information. if (!this->schedule()[0].rescoup.empty()) { auto master_mode = this->schedule()[0].rescoup().masterMode(); if (master_mode) { this->reservoirCouplingMaster_ = std::make_unique( FlowGenericVanguard::comm(), this->schedule() ); this->reservoirCouplingMaster_->spawnSlaveProcesses(argc, argv); } } } simulator_.setEpisodeIndex(-1); // Create timers and file for writing timing info. solverTimer_ = std::make_unique(); totalTimer_ = std::make_unique(); totalTimer_->start(); // adaptive time stepping bool enableAdaptive = Parameters::Get(); bool enableTUNING = Parameters::Get(); if (enableAdaptive) { const UnitSystem& unitSystem = this->simulator_.vanguard().eclState().getUnits(); const auto& sched_state = schedule()[timer.currentStepNum()]; auto max_next_tstep = sched_state.max_next_tstep(enableTUNING); if (enableTUNING) { adaptiveTimeStepping_ = std::make_unique(max_next_tstep, sched_state.tuning(), unitSystem, terminalOutput_); } else { adaptiveTimeStepping_ = std::make_unique(unitSystem, max_next_tstep, terminalOutput_); } if (isRestart()) { // For restarts the simulator may have gotten some information // about the next timestep size from the OPMEXTRA field adaptiveTimeStepping_->setSuggestedNextStep(simulator_.timeStepSize()); } } } void updateTUNING(const Tuning& tuning) { modelParam_.tolerance_mb_ = tuning.XXXMBE; if (terminalOutput_) { OpmLog::debug(fmt::format("Setting SimulatorFullyImplicitBlackoil mass balance limit (XXXMBE) to {:.2e}", tuning.XXXMBE)); } } bool runStep(SimulatorTimer& timer) { if (schedule().exitStatus().has_value()) { if (terminalOutput_) { OpmLog::info("Stopping simulation since EXIT was triggered by an action keyword."); } report_.success.exit_status = schedule().exitStatus().value(); return false; } if (serializer_.shouldLoad()) { serializer_.loadTimerInfo(timer); } // Report timestep. if (terminalOutput_) { std::ostringstream ss; timer.report(ss); OpmLog::debug(ss.str()); } if (terminalOutput_) { details::outputReportStep(timer); } // write the inital state at the report stage if (timer.initialStep()) { Dune::Timer perfTimer; perfTimer.start(); simulator_.setEpisodeIndex(-1); simulator_.setEpisodeLength(0.0); simulator_.setTimeStepSize(0.0); wellModel_().beginReportStep(timer.currentStepNum()); simulator_.problem().writeOutput(true); report_.success.output_write_time += perfTimer.stop(); } // Run a multiple steps of the solver depending on the time step control. solverTimer_->start(); if (!solver_) { solver_ = createSolver(wellModel_()); } simulator_.startNextEpisode( simulator_.startTime() + schedule().seconds(timer.currentStepNum()), timer.currentStepLength()); simulator_.setEpisodeIndex(timer.currentStepNum()); if (serializer_.shouldLoad()) { wellModel_().prepareDeserialize(serializer_.loadStep() - 1); serializer_.loadState(); simulator_.model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0); } this->solver_->model().beginReportStep(); const bool enableTUNING = Parameters::Get(); // If sub stepping is enabled allow the solver to sub cycle // in case the report steps are too large for the solver to converge // // \Note: The report steps are met in any case // \Note: The sub stepping will require a copy of the state variables if (adaptiveTimeStepping_) { auto tuningUpdater = [enableTUNING, this, reportStep = timer.currentStepNum()](const double curr_time, double dt, const int timeStep) { auto& schedule = this->simulator_.vanguard().schedule(); auto& events = this->schedule()[reportStep].events(); bool result = false; if (events.hasEvent(ScheduleEvents::TUNING_CHANGE)) { // Unset the event to not trigger it again on the next sub step schedule.clear_event(ScheduleEvents::TUNING_CHANGE, reportStep); const auto& sched_state = schedule[reportStep]; const auto& max_next_tstep = sched_state.max_next_tstep(enableTUNING); const auto& tuning = sched_state.tuning(); if (enableTUNING) { adaptiveTimeStepping_->updateTUNING(max_next_tstep, tuning); // \Note: Assumes TUNING is only used with adaptive time-stepping // \Note: Need to update both solver (model) and simulator since solver is re-created each report step. solver_->model().updateTUNING(tuning); this->updateTUNING(tuning); dt = this->adaptiveTimeStepping_->suggestedNextStep(); } else { dt = max_next_tstep; this->adaptiveTimeStepping_->updateNEXTSTEP(max_next_tstep); } result = max_next_tstep > 0; } const auto& wcycle = schedule[reportStep].wcycle.get(); if (wcycle.empty()) { return result; } const auto& wmatcher = schedule.wellMatcher(reportStep); double wcycle_time_step = wcycle.nextTimeStep(curr_time, dt, wmatcher, this->wellModel_().wellOpenTimes(), this->wellModel_().wellCloseTimes(), [timeStep, &wg_events = schedule[reportStep].wellgroup_events()](const std::string& name) { if (timeStep != 0) { return false; } return wg_events.hasEvent(name, ScheduleEvents::REQUEST_OPEN_WELL); }); wcycle_time_step = this->grid().comm().min(wcycle_time_step); if (dt != wcycle_time_step) { this->adaptiveTimeStepping_->updateNEXTSTEP(wcycle_time_step); return true; } return result; }; tuningUpdater(timer.simulationTimeElapsed(), this->adaptiveTimeStepping_->suggestedNextStep(), 0); const auto& events = schedule()[timer.currentStepNum()].events(); bool event = events.hasEvent(ScheduleEvents::NEW_WELL) || events.hasEvent(ScheduleEvents::INJECTION_TYPE_CHANGED) || events.hasEvent(ScheduleEvents::WELL_SWITCHED_INJECTOR_PRODUCER) || events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE) || events.hasEvent(ScheduleEvents::INJECTION_UPDATE) || events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE); auto stepReport = adaptiveTimeStepping_->step(timer, *solver_, event, tuningUpdater); report_ += stepReport; //Pass simulation report to eclwriter for summary output simulator_.problem().setSimulationReport(report_); } else { // solve for complete report step auto stepReport = solver_->step(timer); report_ += stepReport; if (terminalOutput_) { std::ostringstream ss; stepReport.reportStep(ss); OpmLog::info(ss.str()); } } // write simulation state at the report stage Dune::Timer perfTimer; perfTimer.start(); const double nextstep = adaptiveTimeStepping_ ? adaptiveTimeStepping_->suggestedNextStep() : -1.0; simulator_.problem().setNextTimeStepSize(nextstep); simulator_.problem().writeOutput(true); report_.success.output_write_time += perfTimer.stop(); solver_->model().endReportStep(); // take time that was used to solve system for this reportStep solverTimer_->stop(); // update timing. report_.success.solver_time += solverTimer_->secsSinceStart(); if (this->grid().comm().rank() == 0) { // Grab the step convergence reports that are new since last we // were here. const auto& reps = this->solver_->model().stepReports(); convergence_output_.write(reps); } // Increment timer, remember well state. ++timer; if (terminalOutput_) { std::string msg = "Time step took " + std::to_string(solverTimer_->secsSinceStart()) + " seconds; " "total solver time " + std::to_string(report_.success.solver_time) + " seconds."; OpmLog::debug(msg); } serializer_.save(timer); return true; } SimulatorReport finalize() { // make sure all output is written to disk before run is finished { Dune::Timer finalOutputTimer; finalOutputTimer.start(); simulator_.problem().finalizeOutput(); report_.success.output_write_time += finalOutputTimer.stop(); } // Stop timer and create timing report totalTimer_->stop(); report_.success.total_time = totalTimer_->secsSinceStart(); report_.success.converged = true; return report_; } const Grid& grid() const { return simulator_.vanguard().grid(); } template void serializeOp(Serializer& serializer) { serializer(simulator_); serializer(report_); serializer(adaptiveTimeStepping_); } const Model& model() const { return solver_->model(); } protected: //! \brief Load simulator state from hdf5 serializer. void loadState([[maybe_unused]] HDF5Serializer& serializer, [[maybe_unused]] const std::string& groupName) override { #if HAVE_HDF5 serializer.read(*this, groupName, "simulator_data"); #endif } //! \brief Save simulator state using hdf5 serializer. void saveState([[maybe_unused]] HDF5Serializer& serializer, [[maybe_unused]] const std::string& groupName) const override { #if HAVE_HDF5 serializer.write(*this, groupName, "simulator_data"); #endif } //! \brief Returns header data std::array getHeader() const override { std::ostringstream str; Parameters::printValues(str); return {"OPM Flow", moduleVersion(), compileTimestamp(), simulator_.vanguard().caseName(), str.str()}; } //! \brief Returns local-to-global cell mapping. const std::vector& getCellMapping() const override { return simulator_.vanguard().globalCell(); } std::unique_ptr createSolver(WellModel& wellModel) { auto model = std::make_unique(simulator_, modelParam_, wellModel, terminalOutput_); if (this->modelParam_.write_partitions_) { const auto& iocfg = this->eclState().cfg().io(); const auto odir = iocfg.getOutputDir() / std::filesystem::path { "partition" } / iocfg.getBaseName(); if (this->grid().comm().rank() == 0) { create_directories(odir); } this->grid().comm().barrier(); model->writePartitions(odir); this->modelParam_.write_partitions_ = false; } return std::make_unique(solverParam_, std::move(model)); } const EclipseState& eclState() const { return simulator_.vanguard().eclState(); } const Schedule& schedule() const { return simulator_.vanguard().schedule(); } bool isRestart() const { const auto& initconfig = eclState().getInitConfig(); return initconfig.restartRequested(); } WellModel& wellModel_() { return simulator_.problem().wellModel(); } const WellModel& wellModel_() const { return simulator_.problem().wellModel(); } // Data. Simulator& simulator_; ModelParameters modelParam_; SolverParameters solverParam_; std::unique_ptr solver_; // Observed objects. PhaseUsage phaseUsage_; // Misc. data bool terminalOutput_; SimulatorReport report_; std::unique_ptr solverTimer_; std::unique_ptr totalTimer_; std::unique_ptr adaptiveTimeStepping_; SimulatorConvergenceOutput convergence_output_; bool slaveMode_{false}; std::unique_ptr reservoirCouplingMaster_{nullptr}; std::unique_ptr reservoirCouplingSlave_{nullptr}; SimulatorSerializer serializer_; }; } // namespace Opm #endif // OPM_SIMULATOR_FULLY_IMPLICIT_BLACKOIL_HEADER_INCLUDED