28#ifndef OPM_ECL_WRITER_HPP
29#define OPM_ECL_WRITER_HPP
31#include <dune/grid/common/partitionset.hh>
33#include <opm/common/TimingMacros.hpp>
34#include <opm/common/OpmLog/OpmLog.hpp>
35#include <opm/input/eclipse/Schedule/RPTConfig.hpp>
37#include <opm/input/eclipse/Units/UnitSystem.hpp>
38#include <opm/input/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp>
40#include <opm/output/eclipse/Inplace.hpp>
41#include <opm/output/eclipse/RestartValue.hpp>
57#ifdef RESERVOIR_COUPLING_ENABLED
61#include <boost/date_time/posix_time/posix_time.hpp>
114template <
class TypeTag,
class OutputModule>
116 GetPropType<TypeTag, Properties::EquilGrid>,
117 GetPropType<TypeTag, Properties::GridView>,
118 GetPropType<TypeTag, Properties::ElementMapper>,
119 GetPropType<TypeTag, Properties::Scalar>>
129 using Element =
typename GridView::template Codim<0>::Entity;
131 using ElementIterator =
typename GridView::template Codim<0>::Iterator;
134 typedef Dune::MultipleCodimMultipleGeomTypeMapper< GridView > VertexMapper;
136 enum { enableEnergy = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal ||
137 getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal };
138 enum { enableMech = getPropValue<TypeTag, Properties::EnableMech>() };
139 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
140 enum { enableGeochemistry = getPropValue<TypeTag, Properties::EnableGeochemistry>() };
146 OutputModule::registerParameters();
148 Parameters::Register<Parameters::EnableAsyncEclOutput>
149 (
"Write the ECL-formated results in a non-blocking way "
150 "(i.e., using a separate thread).");
151 Parameters::Register<Parameters::EnableEsmry>
152 (
"Write ESMRY file for fast loading of summary data.");
159 :
BaseType(simulator.vanguard().schedule(),
160 simulator.vanguard().eclState(),
161 simulator.vanguard().summaryConfig(),
162 simulator.vanguard().grid(),
163 ((simulator.vanguard().grid().comm().rank() == 0)
164 ? &simulator.vanguard().equilGrid()
166 simulator.vanguard().gridView(),
167 simulator.vanguard().cartesianIndexMapper(),
168 ((simulator.vanguard().grid().comm().rank() == 0)
169 ? &simulator.vanguard().equilCartesianIndexMapper()
171 Parameters::
Get<Parameters::EnableAsyncEclOutput>(),
172 Parameters::
Get<Parameters::EnableEsmry>())
173 , simulator_(simulator)
176 if (this->simulator_.vanguard().grid().comm().size() > 1) {
177 auto smryCfg = (this->simulator_.vanguard().grid().comm().rank() == 0)
178 ? this->
eclIO_->finalSummaryConfig()
181 eclBroadcast(this->simulator_.vanguard().grid().comm(), smryCfg);
183 this->outputModule_ = std::make_unique<OutputModule>
189 this->outputModule_ = std::make_unique<OutputModule>
190 (simulator, this->
eclIO_->finalSummaryConfig(), this->collectOnIORank_);
193 this->rank_ = this->simulator_.vanguard().grid().comm().rank();
195 this->simulator_.vanguard().eclState().computeFipRegionStatistics();
203 return simulator_.vanguard().equilGrid();
212 const int reportStepNum = simulator_.episodeIndex() + 1;
232 if (reportStepNum == 0)
235 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
236 const Scalar totalCpuTime =
237 simulator_.executionTimer().realTimeElapsed() +
238 simulator_.setupTimer().realTimeElapsed() +
239 simulator_.vanguard().setupTime();
241 const auto localWellData = simulator_.problem().wellModel().wellData();
242 const auto localWBP = simulator_.problem().wellModel().wellBlockAveragePressures();
243 const auto localGroupAndNetworkData = simulator_.problem().wellModel()
244 .groupAndNetworkData(reportStepNum);
246 const auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
247 const auto localWellTestState = simulator_.problem().wellModel().wellTestState();
248 this->prepareLocalCellData(isSubStep, reportStepNum);
250 if (this->outputModule_->needInterfaceFluxes(isSubStep)) {
251 this->captureLocalFluxData();
257 std::map<std::pair<std::string,int>,
double> dummy;
259 outputModule_->getBlockData(),
263 localGroupAndNetworkData,
266 this->outputModule_->getInterRegFlows(),
273 if (! iregFlows.readIsConsistent()) {
274 throw std::runtime_error {
275 "Inconsistent inter-region flow "
276 "region set names in parallel"
284 this->simulator_.vanguard().grid().comm());
288 std::map<std::string, double> miscSummaryData;
289 std::map<std::string, std::vector<double>> regionData;
293 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
295 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
303 if (totalCpuTime != 0.0) {
304 miscSummaryData[
"TCPU"] = totalCpuTime;
330 const auto rcGroupRates = this->collectReservoirCouplingGroupRates_();
337 : this->outputModule_->getBlockData();
341 : this->outputModule_->getInterRegFlows();
347 localGroupAndNetworkData,
353 this->outputModule_->initialInplace(),
355 this->summaryState(),
357 rcGroupRates ? &(*rcGroupRates) :
nullptr);
364 const auto& gridView = simulator_.vanguard().gridView();
368 this->outputModule_->
369 allocBuffers(num_interior, 0,
false,
false,
false);
372#pragma omp parallel for
374 for (
int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
375 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, 0);
376 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
378 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
383 outputModule_->calc_initial_inplace(simulator_.gridView().comm());
386 const auto& fip = simulator_.vanguard().eclState().getEclipseConfig().fip();
387 if (fip.output(FIPConfig::OutputField::FIELD) ||
388 fip.output(FIPConfig::OutputField::RESV))
390 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
392 const auto start_time = boost::posix_time::
393 from_time_t(simulator_.vanguard().schedule().getStartTime());
396 this->inplace_ = *this->outputModule_->initialInplace();
398 this->outputModule_->
399 outputFipAndResvLog(this->inplace_, 0, 0.0, start_time,
400 false, simulator_.gridView().comm());
404 outputModule_->outputFipAndResvLogToCSV(0,
false, simulator_.gridView().comm());
421 const auto firstStep = this->initialStep();
425 const auto& rpt = this->
schedule_[simStep].rpt_config();
427 if (rpt.contains(
"WELSPECS") && (rpt.at(
"WELSPECS") > 0)) {
430 this->writeWellspecReport(timer);
438 if (rpt.contains(
"WELLS") && rpt.at(
"WELLS") > 0) {
439 this->writeWellflowReport(timer, simStep, rpt.at(
"WELLS"));
442 this->outputModule_->outputFipAndResvLog(this->inplace_,
447 simulator_.gridView().comm());
452 void writeOutput(data::Solution&& localCellData,
const bool isSubStep,
const bool isForcedFinalOutput)
456 const int reportStepNum = simulator_.episodeIndex() + 1;
457 this->prepareLocalCellData(isSubStep, reportStepNum);
458 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
461 auto localWellData = simulator_.problem().wellModel().wellData();
462 auto localGroupAndNetworkData = simulator_.problem().wellModel()
463 .groupAndNetworkData(reportStepNum);
465 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
466 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
468 const bool isFlowsn = this->outputModule_->getFlows().hasFlowsn();
469 auto flowsn = this->outputModule_->getFlows().getFlowsn();
471 const bool isFloresn = this->outputModule_->getFlows().hasFloresn();
472 auto floresn = this->outputModule_->getFlows().getFloresn();
474 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
476 if (localCellData.empty()) {
477 this->outputModule_->assignToSolution(localCellData);
481 this->outputModule_->addRftDataToWells(localWellData,
483 simulator_.gridView().comm());
487 this->collectOnIORank_.doesNeedReordering())
495 this->outputModule_->getBlockData(),
496 this->outputModule_->getExtraBlockData(),
499 localGroupAndNetworkData,
509 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
513 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
514 const Scalar nextStepSize = simulator_.problem().nextTimeStepSize();
515 std::optional<int> timeStepIdx;
516 if (Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
517 timeStepIdx = simulator_.timeStepIndex();
521 std::move(localCellData),
522 std::move(localWellData),
523 std::move(localGroupAndNetworkData),
524 std::move(localAquiferData),
525 std::move(localWellTestState),
528 this->summaryState(),
529 this->simulator_.problem().thresholdPressure().getRestartVector(),
530 curTime, nextStepSize,
531 Parameters::Get<Parameters::EclOutputDoublePrecision>(),
532 isFlowsn, std::move(flowsn),
533 isFloresn, std::move(floresn));
539 const auto enablePCHysteresis = simulator_.problem().materialLawManager()->enablePCHysteresis();
540 const auto enableNonWettingHysteresis = simulator_.problem().materialLawManager()->enableNonWettingHysteresis();
541 const auto enableWettingHysteresis = simulator_.problem().materialLawManager()->enableWettingHysteresis();
542 const auto oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
543 const auto gasActive = FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
544 const auto waterActive = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
545 const auto enableSwatinit = simulator_.vanguard().eclState().fieldProps().has_double(
"SWATINIT");
547 std::vector<RestartKey> solutionKeys {
548 {
"PRESSURE", UnitSystem::measure::pressure},
549 {
"SWAT", UnitSystem::measure::identity, waterActive},
550 {
"SGAS", UnitSystem::measure::identity, gasActive},
551 {
"TEMP", UnitSystem::measure::temperature, enableEnergy},
552 {
"SSOLVENT", UnitSystem::measure::identity, enableSolvent},
554 {
"RS", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGas()},
555 {
"RV", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedOil()},
556 {
"RVW", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedWater()},
557 {
"RSW", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGasInWater()},
559 {
"SGMAX", UnitSystem::measure::identity, enableNonWettingHysteresis && oilActive && gasActive},
560 {
"SHMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && gasActive},
562 {
"SOMAX", UnitSystem::measure::identity,
563 (enableNonWettingHysteresis && oilActive && waterActive)
564 || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
566 {
"SOMIN", UnitSystem::measure::identity, enablePCHysteresis && oilActive && gasActive},
567 {
"SWHY1", UnitSystem::measure::identity, enablePCHysteresis && oilActive && waterActive},
568 {
"SWMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && waterActive},
570 {
"PPCW", UnitSystem::measure::pressure, enableSwatinit},
574 const auto& tracers = simulator_.vanguard().eclState().tracer();
576 for (
const auto& tracer : tracers) {
577 const auto enableSolTracer =
578 ((tracer.phase == Phase::GAS) && FluidSystem::enableDissolvedGas()) ||
579 ((tracer.phase == Phase::OIL) && FluidSystem::enableVaporizedOil());
581 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity,
true);
582 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
586 const auto& inputThpres = eclState().getSimulationConfig().getThresholdPressure();
587 const std::vector<RestartKey> extraKeys {
588 {
"OPMEXTRA", UnitSystem::measure::identity,
false},
589 {
"THRESHPR", UnitSystem::measure::pressure, inputThpres.active()},
592 const auto& gridView = this->simulator_.vanguard().gridView();
593 const auto numElements = gridView.size(0);
597 this->outputModule_->allocBuffers(numElements,
603 const auto restartSolution =
605 solutionKeys, gridView.comm(), 0);
607 if (!restartSolution.empty()) {
608 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
610 this->outputModule_->setRestart(restartSolution, elemIdx, globalIdx);
613 this->simulator_.problem().readSolutionFromOutputModule(0,
true);
614 this->simulator_.problem().temperatureModel().init();
615 ElementContext elemCtx(this->simulator_);
616 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
617 elemCtx.updatePrimaryStencil(elem);
618 elemCtx.updatePrimaryIntensiveQuantities(0);
620 this->outputModule_->updateFluidInPlace(elemCtx);
623 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
631 const auto restartStepIdx = this->simulator_.vanguard()
632 .eclState().getInitConfig().getRestartStep();
634 this->outputModule_->allocBuffers(numElements,
642 const auto restartValues =
645 this->summaryState(),
646 solutionKeys, extraKeys, gridView.comm());
648 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
650 this->outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
653 auto& tracer_model = simulator_.problem().tracerModel();
654 for (
int tracer_index = 0; tracer_index < tracer_model.numTracers(); ++tracer_index) {
657 const auto& free_tracer_name = tracer_model.fname(tracer_index);
658 const auto& free_tracer_solution = restartValues.solution
659 .template data<double>(free_tracer_name);
661 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
663 tracer_model.setFreeTracerConcentration
664 (tracer_index, elemIdx, free_tracer_solution[globalIdx]);
669 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
670 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil()))
672 tracer_model.setEnableSolTracers(tracer_index,
true);
674 const auto& sol_tracer_name = tracer_model.sname(tracer_index);
675 const auto& sol_tracer_solution = restartValues.solution
676 .template data<double>(sol_tracer_name);
678 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
680 tracer_model.setSolTracerConcentration
681 (tracer_index, elemIdx, sol_tracer_solution[globalIdx]);
685 tracer_model.setEnableSolTracers(tracer_index,
false);
687 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
688 tracer_model.setSolTracerConcentration(tracer_index, elemIdx, 0.0);
693 if (inputThpres.active()) {
694 const_cast<Simulator&
>(this->simulator_)
695 .problem().thresholdPressure()
696 .setFromRestart(restartValues.getExtra(
"THRESHPR"));
699 restartTimeStepSize_ = restartValues.getExtra(
"OPMEXTRA")[0];
700 if (restartTimeStepSize_ <= 0) {
701 restartTimeStepSize_ = std::numeric_limits<double>::max();
705 this->simulator_.problem().wellModel()
706 .initFromRestartFile(restartValues);
708 if (!restartValues.aquifer.empty()) {
709 this->simulator_.problem().mutableAquiferModel()
710 .initFromRestart(restartValues.aquifer);
720 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
723 if (
const auto* iip = this->outputModule_->initialInplace(); iip !=
nullptr) {
724 this->inplace_ = *iip;
730 {
return *outputModule_; }
733 {
return *outputModule_; }
736 {
return restartTimeStepSize_; }
738 template <
class Serializer>
741 serializer(*outputModule_);
745 static bool enableEclOutput_()
747 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
751 const EclipseState& eclState()
const
752 {
return simulator_.vanguard().eclState(); }
754 SummaryState& summaryState()
755 {
return simulator_.vanguard().summaryState(); }
757 Action::State& actionState()
758 {
return simulator_.vanguard().actionState(); }
761 {
return simulator_.vanguard().udqState(); }
763 const Schedule& schedule()
const
764 {
return simulator_.vanguard().schedule(); }
768 std::optional<data::ReservoirCouplingGroupRates> collectReservoirCouplingGroupRates_()
770#ifdef RESERVOIR_COUPLING_ENABLED
775 using WellModelType = std::remove_cvref_t<
776 decltype(simulator_.problem().wellModel())>;
777 if constexpr (
requires(WellModelType& wm) { wm.isReservoirCouplingMaster(); }) {
778 auto& wellModel = simulator_.problem().wellModel();
779 if (!wellModel.isReservoirCouplingMaster()) {
782 return wellModel.reservoirCouplingMaster()
783 .collectGroupRatesForSummary();
789 void prepareLocalCellData(
const bool isSubStep,
790 const int reportStepNum)
792 OPM_TIMEBLOCK(prepareLocalCellData);
794 if (this->outputModule_->localDataValid()) {
798 const auto& gridView = simulator_.vanguard().gridView();
803 this->outputModule_->
804 allocBuffers(num_interior, reportStepNum,
805 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
808 ElementContext elemCtx(simulator_);
813 OPM_TIMEBLOCK(prepareCellBasedData);
815 this->outputModule_->prepareDensityAccumulation();
816 this->outputModule_->setupExtractors(isSubStep, reportStepNum);
817 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
818 elemCtx.updatePrimaryStencil(elem);
819 elemCtx.updatePrimaryIntensiveQuantities(0);
821 this->outputModule_->processElement(elemCtx);
822 this->outputModule_->processElementBlockData(elemCtx);
824 this->outputModule_->clearExtractors();
826 this->outputModule_->accumulateDensityParallel();
830 OPM_TIMEBLOCK(prepareFluidInPlace);
833#pragma omp parallel for
835 for (
int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
836 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, 0);
837 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
839 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
843 this->outputModule_->validateLocalData();
846 this->simulator_.vanguard().grid().comm());
849 void captureLocalFluxData()
851 OPM_TIMEBLOCK(captureLocalData);
853 const auto& gridView = this->simulator_.vanguard().gridView();
854 const auto timeIdx = 0u;
856 auto elemCtx = ElementContext { this->simulator_ };
858 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
859 const auto activeIndex = [&elemMapper](
const Element& e)
861 return elemMapper.index(e);
864 const auto cartesianIndex = [
this](
const int elemIndex)
866 return this->
cartMapper_.cartesianIndex(elemIndex);
869 this->outputModule_->initializeFluxData();
873 for (
const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
874 elemCtx.updateStencil(elem);
875 elemCtx.updateIntensiveQuantities(timeIdx);
876 elemCtx.updateExtensiveQuantities(timeIdx);
878 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
882 this->simulator_.vanguard().grid().comm())
884 this->outputModule_->finalizeFluxData();
887 void writeWellspecReport(const SimulatorTimer& timer)
const
889 const auto changedWells = this->
schedule_
890 .changed_wells(timer.reportStepNum(), this->initialStep());
892 const auto changedWellLists = this->
schedule_
893 .changedWellLists(timer.reportStepNum(), this->initialStep());
895 if (changedWells.empty() && !changedWellLists) {
899 this->outputModule_->outputWellspecReport(changedWells,
901 timer.reportStepNum(),
902 timer.simulationTimeElapsed(),
903 timer.currentDateTime());
906 void writeWellflowReport(
const SimulatorTimer& timer,
908 const int wellsRequest)
const
910 this->outputModule_->outputTimeStamp(
"WELLS",
911 timer.simulationTimeElapsed(),
912 timer.reportStepNum(),
913 timer.currentDateTime());
915 const auto wantConnData = wellsRequest > 1;
917 this->outputModule_->outputProdLog(simStep, wantConnData);
918 this->outputModule_->outputInjLog(simStep, wantConnData);
919 this->outputModule_->outputCumLog(simStep, wantConnData);
920 this->outputModule_->outputMSWLog(simStep);
923 int initialStep()
const
925 const auto& initConfig = this->eclState().cfg().init();
927 return initConfig.restartRequested()
928 ? initConfig.getRestartStep()
932 Simulator& simulator_;
933 std::unique_ptr<OutputModule> outputModule_;
934 Scalar restartTimeStepSize_;
#define OPM_END_PARALLEL_TRY_CATCH(prefix, comm)
Catch exception and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:192
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
Contains the classes required to extend the black-oil model by energy.
Declares the properties required by the black oil model.
int localIdxToGlobalIdx(unsigned localIdx) const
Definition: CollectDataOnIORank_impl.hpp:1097
InterRegFlowMap & globalInterRegFlows()
Definition: CollectDataOnIORank.hpp:114
bool isParallel() const
Definition: CollectDataOnIORank.hpp:129
bool isIORank() const
Definition: CollectDataOnIORank.hpp:126
const std::map< std::pair< std::string, int >, double > & globalBlockData() const
Definition: CollectDataOnIORank.hpp:90
void collect(const data::Solution &localCellData, const std::map< std::pair< std::string, int >, double > &localBlockData, std::map< std::pair< std::string, int >, double > &localExtraBlockData, const data::Wells &localWellData, const data::WellBlockAveragePressures &localWBPData, const data::GroupAndNetworkValues &localGroupAndNetworkData, const data::Aquifers &localAquiferData, const WellTestState &localWellTestState, const InterRegFlowMap &interRegFlows, const std::array< FlowsData< double >, 3 > &localFlowsn, const std::array< FlowsData< double >, 3 > &localFloresn)
Definition: CollectDataOnIORank_impl.hpp:964
const data::Solution & globalCellData() const
Definition: CollectDataOnIORank.hpp:93
Definition: EclGenericWriter.hpp:69
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:158
const Schedule & schedule_
Definition: EclGenericWriter.hpp:161
SimulatorReport simulation_report_
Definition: EclGenericWriter.hpp:171
SimulatorReportSingle sub_step_report_
Definition: EclGenericWriter.hpp:170
void evalSummary(int reportStepNum, GetPropType< TypeTag, Properties::Scalar > curTime, const data::Wells &localWellData, const data::WellBlockAveragePressures &localWBPData, const data::GroupAndNetworkValues &localGroupAndNetworkData, const std::map< int, data::AquiferData > &localAquiferData, const std::map< std::pair< std::string, int >, double > &blockData, const std::map< std::string, double > &miscSummaryData, const std::map< std::string, std::vector< double > > ®ionData, const Inplace &inplace, const Inplace *initialInPlace, const InterRegFlowMap &interRegFlows, SummaryState &summaryState, UDQState &udqState, const data::ReservoirCouplingGroupRates *rcGroupRates=nullptr)
Definition: EclGenericWriter_impl.hpp:1021
const Dune::CartesianIndexMapper< GetPropType< TypeTag, Properties::Grid > > & cartMapper_
Definition: EclGenericWriter.hpp:167
std::unique_ptr< EclipseIO > eclIO_
Definition: EclGenericWriter.hpp:163
void doWriteOutput(const int reportStepNum, const std::optional< int > timeStepNum, const bool isSubStep, const bool forcedSimulationFinished, data::Solution &&localCellData, data::Wells &&localWellData, data::GroupAndNetworkValues &&localGroupAndNetworkData, data::Aquifers &&localAquiferData, WellTestState &&localWTestState, const Action::State &actionState, const UDQState &udqState, const SummaryState &summaryState, const std::vector< GetPropType< TypeTag, Properties::Scalar > > &thresholdPressure, GetPropType< TypeTag, Properties::Scalar > curTime, GetPropType< TypeTag, Properties::Scalar > nextStepSize, bool doublePrecision, bool isFlowsn, std::array< FlowsData< double >, 3 > &&flowsn, bool isFloresn, std::array< FlowsData< double >, 3 > &&floresn)
Definition: EclGenericWriter_impl.hpp:911
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:120
OutputModule & mutableOutputModule() const
Definition: EclWriter.hpp:732
const OutputModule & outputModule() const
Definition: EclWriter.hpp:729
void writeOutput(data::Solution &&localCellData, const bool isSubStep, const bool isForcedFinalOutput)
Definition: EclWriter.hpp:452
void evalSummaryState(bool isSubStep)
collect and pass data and pass it to eclIO writer
Definition: EclWriter.hpp:209
static void registerParameters()
Definition: EclWriter.hpp:144
void serializeOp(Serializer &serializer)
Definition: EclWriter.hpp:739
void writeInitialFIPReport()
Writes the initial FIP report as configured in RPTSOL.
Definition: EclWriter.hpp:362
void beginRestart()
Definition: EclWriter.hpp:537
EclWriter(Simulator &simulator)
Definition: EclWriter.hpp:158
void writeReports(const SimulatorTimer &timer)
Definition: EclWriter.hpp:407
void endRestart()
Definition: EclWriter.hpp:715
Scalar restartTimeStepSize() const
Definition: EclWriter.hpp:735
~EclWriter()
Definition: EclWriter.hpp:198
const EquilGrid & globalGrid() const
Definition: EclWriter.hpp:201
virtual int reportStepNum() const
Current report step number. This might differ from currentStepNum in case of sub stepping.
Definition: SimulatorTimerInterface.hpp:109
Definition: SimulatorTimer.hpp:39
virtual boost::posix_time::ptime currentDateTime() const
Return the current time as a posix time object.
double simulationTimeElapsed() const override
Defines the common properties required by the porous medium multi-phase models.
Definition: ActionHandler.hpp:34
Definition: blackoilnewtonmethodparams.hpp:31
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:187
std::size_t countLocalInteriorCellsGridView(const GridView &gridView)
Get the number of local interior cells in a grid view.
Definition: countGlobalCells.hpp:45
Definition: blackoilbioeffectsmodules.hh:45
data::Solution loadParallelRestartSolution(const EclipseIO *eclIO, const std::vector< RestartKey > &solutionKeys, Parallel::Communication comm, const int step)
void eclBroadcast(Parallel::Communication, T &)
RestartValue loadParallelRestart(const EclipseIO *eclIO, Action::State &actionState, SummaryState &summaryState, const std::vector< RestartKey > &solutionKeys, const std::vector< RestartKey > &extraKeys, Parallel::Communication comm)
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
Definition: EclWriter.hpp:78
static constexpr bool value
Definition: EclWriter.hpp:78
Definition: EclWriter.hpp:75
static constexpr bool value
Definition: EclWriter.hpp:75
Definition: EclWriter.hpp:85
static constexpr bool value
Definition: EclWriter.hpp:85
Definition: EclWriter.hpp:82
static constexpr bool value
Definition: EclWriter.hpp:82
SimulatorReportSingle success
Definition: SimulatorReport.hpp:123
unsigned int min_linear_iterations
Definition: SimulatorReport.hpp:52
unsigned int total_newton_iterations
Definition: SimulatorReport.hpp:50
unsigned int max_linear_iterations
Definition: SimulatorReport.hpp:53
unsigned int total_linear_iterations
Definition: SimulatorReport.hpp:51