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>
56#include <boost/date_time/posix_time/posix_time.hpp>
109template <
class TypeTag,
class OutputModule>
111 GetPropType<TypeTag, Properties::EquilGrid>,
112 GetPropType<TypeTag, Properties::GridView>,
113 GetPropType<TypeTag, Properties::ElementMapper>,
114 GetPropType<TypeTag, Properties::Scalar>>
124 using Element =
typename GridView::template Codim<0>::Entity;
126 using ElementIterator =
typename GridView::template Codim<0>::Iterator;
129 typedef Dune::MultipleCodimMultipleGeomTypeMapper< GridView > VertexMapper;
131 enum { enableEnergy = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal ||
132 getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal };
133 enum { enableMech = getPropValue<TypeTag, Properties::EnableMech>() };
134 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
135 enum { enableGeochemistry = getPropValue<TypeTag, Properties::EnableGeochemistry>() };
141 OutputModule::registerParameters();
143 Parameters::Register<Parameters::EnableAsyncEclOutput>
144 (
"Write the ECL-formated results in a non-blocking way "
145 "(i.e., using a separate thread).");
146 Parameters::Register<Parameters::EnableEsmry>
147 (
"Write ESMRY file for fast loading of summary data.");
154 :
BaseType(simulator.vanguard().schedule(),
155 simulator.vanguard().eclState(),
156 simulator.vanguard().summaryConfig(),
157 simulator.vanguard().grid(),
158 ((simulator.vanguard().grid().comm().rank() == 0)
159 ? &simulator.vanguard().equilGrid()
161 simulator.vanguard().gridView(),
162 simulator.vanguard().cartesianIndexMapper(),
163 ((simulator.vanguard().grid().comm().rank() == 0)
164 ? &simulator.vanguard().equilCartesianIndexMapper()
166 Parameters::
Get<Parameters::EnableAsyncEclOutput>(),
167 Parameters::
Get<Parameters::EnableEsmry>())
168 , simulator_(simulator)
171 if (this->simulator_.vanguard().grid().comm().size() > 1) {
172 auto smryCfg = (this->simulator_.vanguard().grid().comm().rank() == 0)
173 ? this->
eclIO_->finalSummaryConfig()
176 eclBroadcast(this->simulator_.vanguard().grid().comm(), smryCfg);
178 this->outputModule_ = std::make_unique<OutputModule>
184 this->outputModule_ = std::make_unique<OutputModule>
185 (simulator, this->
eclIO_->finalSummaryConfig(), this->collectOnIORank_);
188 this->rank_ = this->simulator_.vanguard().grid().comm().rank();
190 this->simulator_.vanguard().eclState().computeFipRegionStatistics();
198 return simulator_.vanguard().equilGrid();
207 const int reportStepNum = simulator_.episodeIndex() + 1;
227 if (reportStepNum == 0)
230 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
231 const Scalar totalCpuTime =
232 simulator_.executionTimer().realTimeElapsed() +
233 simulator_.setupTimer().realTimeElapsed() +
234 simulator_.vanguard().setupTime();
236 const auto localWellData = simulator_.problem().wellModel().wellData();
237 const auto localWBP = simulator_.problem().wellModel().wellBlockAveragePressures();
238 const auto localGroupAndNetworkData = simulator_.problem().wellModel()
239 .groupAndNetworkData(reportStepNum);
241 const auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
242 const auto localWellTestState = simulator_.problem().wellModel().wellTestState();
243 this->prepareLocalCellData(isSubStep, reportStepNum);
245 if (this->outputModule_->needInterfaceFluxes(isSubStep)) {
246 this->captureLocalFluxData();
252 std::map<std::pair<std::string,int>,
double> dummy;
254 outputModule_->getBlockData(),
258 localGroupAndNetworkData,
261 this->outputModule_->getInterRegFlows(),
268 if (! iregFlows.readIsConsistent()) {
269 throw std::runtime_error {
270 "Inconsistent inter-region flow "
271 "region set names in parallel"
279 this->simulator_.vanguard().grid().comm());
283 std::map<std::string, double> miscSummaryData;
284 std::map<std::string, std::vector<double>> regionData;
288 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
290 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
298 if (totalCpuTime != 0.0) {
299 miscSummaryData[
"TCPU"] = totalCpuTime;
328 : this->outputModule_->getBlockData();
332 : this->outputModule_->getInterRegFlows();
338 localGroupAndNetworkData,
344 this->outputModule_->initialInplace(),
346 this->summaryState(),
354 const auto& gridView = simulator_.vanguard().gridView();
358 this->outputModule_->
359 allocBuffers(num_interior, 0,
false,
false,
false);
362#pragma omp parallel for
364 for (
int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
365 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, 0);
366 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
368 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
373 outputModule_->calc_initial_inplace(simulator_.gridView().comm());
376 const auto& fip = simulator_.vanguard().eclState().getEclipseConfig().fip();
377 if (fip.output(FIPConfig::OutputField::FIELD) ||
378 fip.output(FIPConfig::OutputField::RESV))
380 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
382 const auto start_time = boost::posix_time::
383 from_time_t(simulator_.vanguard().schedule().getStartTime());
386 this->inplace_ = *this->outputModule_->initialInplace();
388 this->outputModule_->
389 outputFipAndResvLog(this->inplace_, 0, 0.0, start_time,
390 false, simulator_.gridView().comm());
394 outputModule_->outputFipAndResvLogToCSV(0,
false, simulator_.gridView().comm());
411 const auto firstStep = this->initialStep();
415 const auto& rpt = this->
schedule_[simStep].rpt_config();
417 if (rpt.contains(
"WELSPECS") && (rpt.at(
"WELSPECS") > 0)) {
420 this->writeWellspecReport(timer);
428 if (rpt.contains(
"WELLS") && rpt.at(
"WELLS") > 0) {
429 this->writeWellflowReport(timer, simStep, rpt.at(
"WELLS"));
432 this->outputModule_->outputFipAndResvLog(this->inplace_,
437 simulator_.gridView().comm());
446 const int reportStepNum = simulator_.episodeIndex() + 1;
447 this->prepareLocalCellData(isSubStep, reportStepNum);
448 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
451 auto localWellData = simulator_.problem().wellModel().wellData();
452 auto localGroupAndNetworkData = simulator_.problem().wellModel()
453 .groupAndNetworkData(reportStepNum);
455 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
456 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
458 const bool isFlowsn = this->outputModule_->getFlows().hasFlowsn();
459 auto flowsn = this->outputModule_->getFlows().getFlowsn();
461 const bool isFloresn = this->outputModule_->getFlows().hasFloresn();
462 auto floresn = this->outputModule_->getFlows().getFloresn();
464 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
466 if (localCellData.empty()) {
467 this->outputModule_->assignToSolution(localCellData);
471 this->outputModule_->addRftDataToWells(localWellData,
473 simulator_.gridView().comm());
477 this->collectOnIORank_.doesNeedReordering())
485 this->outputModule_->getBlockData(),
486 this->outputModule_->getExtraBlockData(),
489 localGroupAndNetworkData,
499 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
503 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
504 const Scalar nextStepSize = simulator_.problem().nextTimeStepSize();
505 std::optional<int> timeStepIdx;
506 if (Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
507 timeStepIdx = simulator_.timeStepIndex();
510 std::move(localCellData),
511 std::move(localWellData),
512 std::move(localGroupAndNetworkData),
513 std::move(localAquiferData),
514 std::move(localWellTestState),
517 this->summaryState(),
518 this->simulator_.problem().thresholdPressure().getRestartVector(),
519 curTime, nextStepSize,
520 Parameters::Get<Parameters::EclOutputDoublePrecision>(),
521 isFlowsn, std::move(flowsn),
522 isFloresn, std::move(floresn));
528 const auto enablePCHysteresis = simulator_.problem().materialLawManager()->enablePCHysteresis();
529 const auto enableNonWettingHysteresis = simulator_.problem().materialLawManager()->enableNonWettingHysteresis();
530 const auto enableWettingHysteresis = simulator_.problem().materialLawManager()->enableWettingHysteresis();
531 const auto oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
532 const auto gasActive = FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
533 const auto waterActive = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
534 const auto enableSwatinit = simulator_.vanguard().eclState().fieldProps().has_double(
"SWATINIT");
536 std::vector<RestartKey> solutionKeys {
537 {
"PRESSURE", UnitSystem::measure::pressure},
538 {
"SWAT", UnitSystem::measure::identity, waterActive},
539 {
"SGAS", UnitSystem::measure::identity, gasActive},
540 {
"TEMP", UnitSystem::measure::temperature, enableEnergy},
541 {
"SSOLVENT", UnitSystem::measure::identity, enableSolvent},
543 {
"RS", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGas()},
544 {
"RV", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedOil()},
545 {
"RVW", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedWater()},
546 {
"RSW", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGasInWater()},
548 {
"SGMAX", UnitSystem::measure::identity, enableNonWettingHysteresis && oilActive && gasActive},
549 {
"SHMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && gasActive},
551 {
"SOMAX", UnitSystem::measure::identity,
552 (enableNonWettingHysteresis && oilActive && waterActive)
553 || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
555 {
"SOMIN", UnitSystem::measure::identity, enablePCHysteresis && oilActive && gasActive},
556 {
"SWHY1", UnitSystem::measure::identity, enablePCHysteresis && oilActive && waterActive},
557 {
"SWMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && waterActive},
559 {
"PPCW", UnitSystem::measure::pressure, enableSwatinit},
563 const auto& tracers = simulator_.vanguard().eclState().tracer();
565 for (
const auto& tracer : tracers) {
566 const auto enableSolTracer =
567 ((tracer.phase == Phase::GAS) && FluidSystem::enableDissolvedGas()) ||
568 ((tracer.phase == Phase::OIL) && FluidSystem::enableVaporizedOil());
570 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity,
true);
571 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
575 const auto& inputThpres = eclState().getSimulationConfig().getThresholdPressure();
576 const std::vector<RestartKey> extraKeys {
577 {
"OPMEXTRA", UnitSystem::measure::identity,
false},
578 {
"THRESHPR", UnitSystem::measure::pressure, inputThpres.active()},
581 const auto& gridView = this->simulator_.vanguard().gridView();
582 const auto numElements = gridView.size(0);
586 this->outputModule_->allocBuffers(numElements,
592 const auto restartSolution =
594 solutionKeys, gridView.comm(), 0);
596 if (!restartSolution.empty()) {
597 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
599 this->outputModule_->setRestart(restartSolution, elemIdx, globalIdx);
602 this->simulator_.problem().readSolutionFromOutputModule(0,
true);
603 this->simulator_.problem().temperatureModel().init();
604 ElementContext elemCtx(this->simulator_);
605 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
606 elemCtx.updatePrimaryStencil(elem);
607 elemCtx.updatePrimaryIntensiveQuantities(0);
609 this->outputModule_->updateFluidInPlace(elemCtx);
612 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
620 const auto restartStepIdx = this->simulator_.vanguard()
621 .eclState().getInitConfig().getRestartStep();
623 this->outputModule_->allocBuffers(numElements,
631 const auto restartValues =
634 this->summaryState(),
635 solutionKeys, extraKeys, gridView.comm());
637 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
639 this->outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
642 auto& tracer_model = simulator_.problem().tracerModel();
643 for (
int tracer_index = 0; tracer_index < tracer_model.numTracers(); ++tracer_index) {
646 const auto& free_tracer_name = tracer_model.fname(tracer_index);
647 const auto& free_tracer_solution = restartValues.solution
648 .template data<double>(free_tracer_name);
650 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
652 tracer_model.setFreeTracerConcentration
653 (tracer_index, elemIdx, free_tracer_solution[globalIdx]);
658 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
659 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil()))
661 tracer_model.setEnableSolTracers(tracer_index,
true);
663 const auto& sol_tracer_name = tracer_model.sname(tracer_index);
664 const auto& sol_tracer_solution = restartValues.solution
665 .template data<double>(sol_tracer_name);
667 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
669 tracer_model.setSolTracerConcentration
670 (tracer_index, elemIdx, sol_tracer_solution[globalIdx]);
674 tracer_model.setEnableSolTracers(tracer_index,
false);
676 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
677 tracer_model.setSolTracerConcentration(tracer_index, elemIdx, 0.0);
682 if (inputThpres.active()) {
683 const_cast<Simulator&
>(this->simulator_)
684 .problem().thresholdPressure()
685 .setFromRestart(restartValues.getExtra(
"THRESHPR"));
688 restartTimeStepSize_ = restartValues.getExtra(
"OPMEXTRA")[0];
689 if (restartTimeStepSize_ <= 0) {
690 restartTimeStepSize_ = std::numeric_limits<double>::max();
694 this->simulator_.problem().wellModel()
695 .initFromRestartFile(restartValues);
697 if (!restartValues.aquifer.empty()) {
698 this->simulator_.problem().mutableAquiferModel()
699 .initFromRestart(restartValues.aquifer);
709 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
712 if (
const auto* iip = this->outputModule_->initialInplace(); iip !=
nullptr) {
713 this->inplace_ = *iip;
719 {
return *outputModule_; }
722 {
return *outputModule_; }
725 {
return restartTimeStepSize_; }
727 template <
class Serializer>
730 serializer(*outputModule_);
734 static bool enableEclOutput_()
736 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
740 const EclipseState& eclState()
const
741 {
return simulator_.vanguard().eclState(); }
743 SummaryState& summaryState()
744 {
return simulator_.vanguard().summaryState(); }
746 Action::State& actionState()
747 {
return simulator_.vanguard().actionState(); }
750 {
return simulator_.vanguard().udqState(); }
752 const Schedule& schedule()
const
753 {
return simulator_.vanguard().schedule(); }
755 void prepareLocalCellData(
const bool isSubStep,
756 const int reportStepNum)
758 OPM_TIMEBLOCK(prepareLocalCellData);
760 if (this->outputModule_->localDataValid()) {
764 const auto& gridView = simulator_.vanguard().gridView();
769 this->outputModule_->
770 allocBuffers(num_interior, reportStepNum,
771 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
774 ElementContext elemCtx(simulator_);
779 OPM_TIMEBLOCK(prepareCellBasedData);
781 this->outputModule_->prepareDensityAccumulation();
782 this->outputModule_->setupExtractors(isSubStep, reportStepNum);
783 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
784 elemCtx.updatePrimaryStencil(elem);
785 elemCtx.updatePrimaryIntensiveQuantities(0);
787 this->outputModule_->processElement(elemCtx);
788 this->outputModule_->processElementBlockData(elemCtx);
790 this->outputModule_->clearExtractors();
792 this->outputModule_->accumulateDensityParallel();
796 OPM_TIMEBLOCK(prepareFluidInPlace);
799#pragma omp parallel for
801 for (
int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
802 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, 0);
803 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
805 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
809 this->outputModule_->validateLocalData();
812 this->simulator_.vanguard().grid().comm());
815 void captureLocalFluxData()
817 OPM_TIMEBLOCK(captureLocalData);
819 const auto& gridView = this->simulator_.vanguard().gridView();
820 const auto timeIdx = 0u;
822 auto elemCtx = ElementContext { this->simulator_ };
824 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
825 const auto activeIndex = [&elemMapper](
const Element& e)
827 return elemMapper.index(e);
830 const auto cartesianIndex = [
this](
const int elemIndex)
832 return this->
cartMapper_.cartesianIndex(elemIndex);
835 this->outputModule_->initializeFluxData();
839 for (
const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
840 elemCtx.updateStencil(elem);
841 elemCtx.updateIntensiveQuantities(timeIdx);
842 elemCtx.updateExtensiveQuantities(timeIdx);
844 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
848 this->simulator_.vanguard().grid().comm())
850 this->outputModule_->finalizeFluxData();
853 void writeWellspecReport(const SimulatorTimer& timer)
const
855 const auto changedWells = this->
schedule_
856 .changed_wells(timer.reportStepNum(), this->initialStep());
858 const auto changedWellLists = this->
schedule_
859 .changedWellLists(timer.reportStepNum(), this->initialStep());
861 if (changedWells.empty() && !changedWellLists) {
865 this->outputModule_->outputWellspecReport(changedWells,
867 timer.reportStepNum(),
868 timer.simulationTimeElapsed(),
869 timer.currentDateTime());
872 void writeWellflowReport(
const SimulatorTimer& timer,
874 const int wellsRequest)
const
876 this->outputModule_->outputTimeStamp(
"WELLS",
877 timer.simulationTimeElapsed(),
878 timer.reportStepNum(),
879 timer.currentDateTime());
881 const auto wantConnData = wellsRequest > 1;
883 this->outputModule_->outputProdLog(simStep, wantConnData);
884 this->outputModule_->outputInjLog(simStep, wantConnData);
885 this->outputModule_->outputCumLog(simStep, wantConnData);
886 this->outputModule_->outputMSWLog(simStep);
889 int initialStep()
const
891 const auto& initConfig = this->eclState().cfg().init();
893 return initConfig.restartRequested()
894 ? initConfig.getRestartStep()
898 Simulator& simulator_;
899 std::unique_ptr<OutputModule> outputModule_;
900 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:67
void doWriteOutput(const int reportStepNum, const std::optional< int > timeStepNum, const bool isSubStep, 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:894
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:154
const Schedule & schedule_
Definition: EclGenericWriter.hpp:157
SimulatorReport simulation_report_
Definition: EclGenericWriter.hpp:167
SimulatorReportSingle sub_step_report_
Definition: EclGenericWriter.hpp:166
const Dune::CartesianIndexMapper< GetPropType< TypeTag, Properties::Grid > > & cartMapper_
Definition: EclGenericWriter.hpp:163
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)
Definition: EclGenericWriter_impl.hpp:1002
std::unique_ptr< EclipseIO > eclIO_
Definition: EclGenericWriter.hpp:159
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:115
OutputModule & mutableOutputModule() const
Definition: EclWriter.hpp:721
const OutputModule & outputModule() const
Definition: EclWriter.hpp:718
void evalSummaryState(bool isSubStep)
collect and pass data and pass it to eclIO writer
Definition: EclWriter.hpp:204
static void registerParameters()
Definition: EclWriter.hpp:139
void serializeOp(Serializer &serializer)
Definition: EclWriter.hpp:728
void writeInitialFIPReport()
Writes the initial FIP report as configured in RPTSOL.
Definition: EclWriter.hpp:352
void beginRestart()
Definition: EclWriter.hpp:526
EclWriter(Simulator &simulator)
Definition: EclWriter.hpp:153
void writeOutput(data::Solution &&localCellData, bool isSubStep)
Definition: EclWriter.hpp:442
void writeReports(const SimulatorTimer &timer)
Definition: EclWriter.hpp:397
void endRestart()
Definition: EclWriter.hpp:704
Scalar restartTimeStepSize() const
Definition: EclWriter.hpp:724
~EclWriter()
Definition: EclWriter.hpp:193
const EquilGrid & globalGrid() const
Definition: EclWriter.hpp:196
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:73
static constexpr bool value
Definition: EclWriter.hpp:73
Definition: EclWriter.hpp:70
static constexpr bool value
Definition: EclWriter.hpp:70
Definition: EclWriter.hpp:80
static constexpr bool value
Definition: EclWriter.hpp:80
Definition: EclWriter.hpp:77
static constexpr bool value
Definition: EclWriter.hpp:77
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