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 enum { enableMech = getPropValue<TypeTag, Properties::EnableMech>() };
133 enum { enableTemperature = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::ConstantTemperature };
134 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
140 OutputModule::registerParameters();
142 Parameters::Register<Parameters::EnableAsyncEclOutput>
143 (
"Write the ECL-formated results in a non-blocking way "
144 "(i.e., using a separate thread).");
145 Parameters::Register<Parameters::EnableEsmry>
146 (
"Write ESMRY file for fast loading of summary data.");
153 :
BaseType(simulator.vanguard().schedule(),
154 simulator.vanguard().eclState(),
155 simulator.vanguard().summaryConfig(),
156 simulator.vanguard().grid(),
157 ((simulator.vanguard().grid().comm().rank() == 0)
158 ? &simulator.vanguard().equilGrid()
160 simulator.vanguard().gridView(),
161 simulator.vanguard().cartesianIndexMapper(),
162 ((simulator.vanguard().grid().comm().rank() == 0)
163 ? &simulator.vanguard().equilCartesianIndexMapper()
165 Parameters::
Get<Parameters::EnableAsyncEclOutput>(),
166 Parameters::
Get<Parameters::EnableEsmry>())
167 , simulator_(simulator)
170 if (this->simulator_.vanguard().grid().comm().size() > 1) {
171 auto smryCfg = (this->simulator_.vanguard().grid().comm().rank() == 0)
172 ? this->
eclIO_->finalSummaryConfig()
175 eclBroadcast(this->simulator_.vanguard().grid().comm(), smryCfg);
177 this->outputModule_ = std::make_unique<OutputModule>
183 this->outputModule_ = std::make_unique<OutputModule>
184 (simulator, this->
eclIO_->finalSummaryConfig(), this->collectOnIORank_);
187 this->rank_ = this->simulator_.vanguard().grid().comm().rank();
189 this->simulator_.vanguard().eclState().computeFipRegionStatistics();
197 return simulator_.vanguard().equilGrid();
206 const int reportStepNum = simulator_.episodeIndex() + 1;
226 if (reportStepNum == 0)
229 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
230 const Scalar totalCpuTime =
231 simulator_.executionTimer().realTimeElapsed() +
232 simulator_.setupTimer().realTimeElapsed() +
233 simulator_.vanguard().setupTime();
235 const auto localWellData = simulator_.problem().wellModel().wellData();
236 const auto localWBP = simulator_.problem().wellModel().wellBlockAveragePressures();
237 const auto localGroupAndNetworkData = simulator_.problem().wellModel()
238 .groupAndNetworkData(reportStepNum);
240 const auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
241 const auto localWellTestState = simulator_.problem().wellModel().wellTestState();
242 this->prepareLocalCellData(isSubStep, reportStepNum);
244 if (this->outputModule_->needInterfaceFluxes(isSubStep)) {
245 this->captureLocalFluxData();
251 std::map<std::pair<std::string,int>,
double> dummy;
253 outputModule_->getBlockData(),
257 localGroupAndNetworkData,
260 this->outputModule_->getInterRegFlows(),
267 if (! iregFlows.readIsConsistent()) {
268 throw std::runtime_error {
269 "Inconsistent inter-region flow "
270 "region set names in parallel"
278 this->simulator_.vanguard().grid().comm());
282 std::map<std::string, double> miscSummaryData;
283 std::map<std::string, std::vector<double>> regionData;
287 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
289 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
297 if (totalCpuTime != 0.0) {
298 miscSummaryData[
"TCPU"] = totalCpuTime;
327 : this->outputModule_->getBlockData();
331 : this->outputModule_->getInterRegFlows();
337 localGroupAndNetworkData,
343 this->outputModule_->initialInplace(),
345 this->summaryState(),
353 const auto& gridView = simulator_.vanguard().gridView();
357 this->outputModule_->
358 allocBuffers(num_interior, 0,
false,
false,
false);
361#pragma omp parallel for
363 for (
int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
364 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, 0);
365 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
367 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
372 outputModule_->calc_initial_inplace(simulator_.gridView().comm());
375 const auto& fip = simulator_.vanguard().eclState().getEclipseConfig().fip();
376 if (fip.output(FIPConfig::OutputField::FIELD) ||
377 fip.output(FIPConfig::OutputField::RESV))
379 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
380 boost::posix_time::ptime start_time =
381 boost::posix_time::from_time_t(simulator_.vanguard().schedule().getStartTime());
384 inplace_ = outputModule_->initialInplace().value();
385 outputModule_->outputFipAndResvLog(inplace_, 0, 0.0, start_time,
386 false, simulator_.gridView().comm());
390 outputModule_->outputFipAndResvLogToCSV(0,
false, simulator_.gridView().comm());
407 const auto firstStep = this->initialStep();
411 const auto& rpt = this->
schedule_[simStep].rpt_config();
413 if (rpt.contains(
"WELSPECS") && (rpt.at(
"WELSPECS") > 0)) {
416 this->writeWellspecReport(timer);
424 if (rpt.contains(
"WELLS") && rpt.at(
"WELLS") > 0) {
425 this->writeWellflowReport(timer, simStep, rpt.at(
"WELLS"));
428 this->outputModule_->outputFipAndResvLog(this->inplace_,
433 simulator_.gridView().comm());
442 const int reportStepNum = simulator_.episodeIndex() + 1;
443 this->prepareLocalCellData(isSubStep, reportStepNum);
444 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
447 auto localWellData = simulator_.problem().wellModel().wellData();
448 auto localGroupAndNetworkData = simulator_.problem().wellModel()
449 .groupAndNetworkData(reportStepNum);
451 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
452 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
454 const bool isFlowsn = this->outputModule_->getFlows().hasFlowsn();
455 auto flowsn = this->outputModule_->getFlows().getFlowsn();
457 const bool isFloresn = this->outputModule_->getFlows().hasFloresn();
458 auto floresn = this->outputModule_->getFlows().getFloresn();
460 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
462 if (localCellData.empty()) {
463 this->outputModule_->assignToSolution(localCellData);
467 this->outputModule_->addRftDataToWells(localWellData,
469 simulator_.gridView().comm());
473 this->collectOnIORank_.doesNeedReordering())
481 this->outputModule_->getBlockData(),
482 this->outputModule_->getExtraBlockData(),
485 localGroupAndNetworkData,
495 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
499 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
500 const Scalar nextStepSize = simulator_.problem().nextTimeStepSize();
501 std::optional<int> timeStepIdx;
502 if (Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
503 timeStepIdx = simulator_.timeStepIndex();
506 std::move(localCellData),
507 std::move(localWellData),
508 std::move(localGroupAndNetworkData),
509 std::move(localAquiferData),
510 std::move(localWellTestState),
513 this->summaryState(),
514 this->simulator_.problem().thresholdPressure().getRestartVector(),
515 curTime, nextStepSize,
516 Parameters::Get<Parameters::EclOutputDoublePrecision>(),
517 isFlowsn, std::move(flowsn),
518 isFloresn, std::move(floresn));
524 const auto enablePCHysteresis = simulator_.problem().materialLawManager()->enablePCHysteresis();
525 const auto enableNonWettingHysteresis = simulator_.problem().materialLawManager()->enableNonWettingHysteresis();
526 const auto enableWettingHysteresis = simulator_.problem().materialLawManager()->enableWettingHysteresis();
527 const auto oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
528 const auto gasActive = FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
529 const auto waterActive = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
530 const auto enableSwatinit = simulator_.vanguard().eclState().fieldProps().has_double(
"SWATINIT");
532 std::vector<RestartKey> solutionKeys {
533 {
"PRESSURE", UnitSystem::measure::pressure},
534 {
"SWAT", UnitSystem::measure::identity, waterActive},
535 {
"SGAS", UnitSystem::measure::identity, gasActive},
536 {
"TEMP", UnitSystem::measure::temperature, enableEnergy},
537 {
"SSOLVENT", UnitSystem::measure::identity, enableSolvent},
539 {
"RS", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGas()},
540 {
"RV", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedOil()},
541 {
"RVW", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedWater()},
542 {
"RSW", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGasInWater()},
544 {
"SGMAX", UnitSystem::measure::identity, enableNonWettingHysteresis && oilActive && gasActive},
545 {
"SHMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && gasActive},
547 {
"SOMAX", UnitSystem::measure::identity,
548 (enableNonWettingHysteresis && oilActive && waterActive)
549 || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
551 {
"SOMIN", UnitSystem::measure::identity, enablePCHysteresis && oilActive && gasActive},
552 {
"SWHY1", UnitSystem::measure::identity, enablePCHysteresis && oilActive && waterActive},
553 {
"SWMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && waterActive},
555 {
"PPCW", UnitSystem::measure::pressure, enableSwatinit},
559 const auto& tracers = simulator_.vanguard().eclState().tracer();
561 for (
const auto& tracer : tracers) {
562 const auto enableSolTracer =
563 ((tracer.phase == Phase::GAS) && FluidSystem::enableDissolvedGas()) ||
564 ((tracer.phase == Phase::OIL) && FluidSystem::enableVaporizedOil());
566 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity,
true);
567 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
571 const auto& inputThpres = eclState().getSimulationConfig().getThresholdPressure();
572 const std::vector<RestartKey> extraKeys {
573 {
"OPMEXTRA", UnitSystem::measure::identity,
false},
574 {
"THRESHPR", UnitSystem::measure::pressure, inputThpres.active()},
577 const auto& gridView = this->simulator_.vanguard().gridView();
578 const auto numElements = gridView.size(0);
582 this->outputModule_->allocBuffers(numElements,
588 const auto restartSolution =
590 solutionKeys, gridView.comm(), 0);
592 if (!restartSolution.empty()) {
593 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
595 this->outputModule_->setRestart(restartSolution, elemIdx, globalIdx);
598 this->simulator_.problem().readSolutionFromOutputModule(0,
true);
599 ElementContext elemCtx(this->simulator_);
600 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
601 elemCtx.updatePrimaryStencil(elem);
602 elemCtx.updatePrimaryIntensiveQuantities(0);
604 this->outputModule_->updateFluidInPlace(elemCtx);
607 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
615 const auto restartStepIdx = this->simulator_.vanguard()
616 .eclState().getInitConfig().getRestartStep();
618 this->outputModule_->allocBuffers(numElements,
626 const auto restartValues =
629 this->summaryState(),
630 solutionKeys, extraKeys, gridView.comm());
632 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
634 this->outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
637 auto& tracer_model = simulator_.problem().tracerModel();
638 for (
int tracer_index = 0; tracer_index < tracer_model.numTracers(); ++tracer_index) {
641 const auto& free_tracer_name = tracer_model.fname(tracer_index);
642 const auto& free_tracer_solution = restartValues.solution
643 .template data<double>(free_tracer_name);
645 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
647 tracer_model.setFreeTracerConcentration
648 (tracer_index, elemIdx, free_tracer_solution[globalIdx]);
653 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
654 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil()))
656 tracer_model.setEnableSolTracers(tracer_index,
true);
658 const auto& sol_tracer_name = tracer_model.sname(tracer_index);
659 const auto& sol_tracer_solution = restartValues.solution
660 .template data<double>(sol_tracer_name);
662 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
664 tracer_model.setSolTracerConcentration
665 (tracer_index, elemIdx, sol_tracer_solution[globalIdx]);
669 tracer_model.setEnableSolTracers(tracer_index,
false);
671 for (
auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
672 tracer_model.setSolTracerConcentration(tracer_index, elemIdx, 0.0);
677 if (inputThpres.active()) {
678 const_cast<Simulator&
>(this->simulator_)
679 .problem().thresholdPressure()
680 .setFromRestart(restartValues.getExtra(
"THRESHPR"));
683 restartTimeStepSize_ = restartValues.getExtra(
"OPMEXTRA")[0];
684 if (restartTimeStepSize_ <= 0) {
685 restartTimeStepSize_ = std::numeric_limits<double>::max();
689 this->simulator_.problem().wellModel()
690 .initFromRestartFile(restartValues);
692 if (!restartValues.aquifer.empty()) {
693 this->simulator_.problem().mutableAquiferModel()
694 .initFromRestart(restartValues.aquifer);
704 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
707 if (this->outputModule_->initialInplace().has_value()) {
708 this->inplace_ = this->outputModule_->initialInplace().value();
714 {
return *outputModule_; }
717 {
return *outputModule_; }
720 {
return restartTimeStepSize_; }
722 template <
class Serializer>
725 serializer(*outputModule_);
729 static bool enableEclOutput_()
731 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
735 const EclipseState& eclState()
const
736 {
return simulator_.vanguard().eclState(); }
738 SummaryState& summaryState()
739 {
return simulator_.vanguard().summaryState(); }
741 Action::State& actionState()
742 {
return simulator_.vanguard().actionState(); }
745 {
return simulator_.vanguard().udqState(); }
747 const Schedule& schedule()
const
748 {
return simulator_.vanguard().schedule(); }
750 void prepareLocalCellData(
const bool isSubStep,
751 const int reportStepNum)
753 OPM_TIMEBLOCK(prepareLocalCellData);
755 if (this->outputModule_->localDataValid()) {
759 const auto& gridView = simulator_.vanguard().gridView();
764 this->outputModule_->
765 allocBuffers(num_interior, reportStepNum,
766 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
769 ElementContext elemCtx(simulator_);
774 OPM_TIMEBLOCK(prepareCellBasedData);
776 this->outputModule_->prepareDensityAccumulation();
777 this->outputModule_->setupExtractors(isSubStep, reportStepNum);
778 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
779 elemCtx.updatePrimaryStencil(elem);
780 elemCtx.updatePrimaryIntensiveQuantities(0);
782 this->outputModule_->processElement(elemCtx);
783 this->outputModule_->processElementBlockData(elemCtx);
785 this->outputModule_->clearExtractors();
787 this->outputModule_->accumulateDensityParallel();
791 OPM_TIMEBLOCK(prepareFluidInPlace);
794#pragma omp parallel for
796 for (
int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
797 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, 0);
798 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
800 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
804 this->outputModule_->validateLocalData();
807 this->simulator_.vanguard().grid().comm());
810 void captureLocalFluxData()
812 OPM_TIMEBLOCK(captureLocalData);
814 const auto& gridView = this->simulator_.vanguard().gridView();
815 const auto timeIdx = 0u;
817 auto elemCtx = ElementContext { this->simulator_ };
819 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
820 const auto activeIndex = [&elemMapper](
const Element& e)
822 return elemMapper.index(e);
825 const auto cartesianIndex = [
this](
const int elemIndex)
827 return this->
cartMapper_.cartesianIndex(elemIndex);
830 this->outputModule_->initializeFluxData();
834 for (
const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
835 elemCtx.updateStencil(elem);
836 elemCtx.updateIntensiveQuantities(timeIdx);
837 elemCtx.updateExtensiveQuantities(timeIdx);
839 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
843 this->simulator_.vanguard().grid().comm())
845 this->outputModule_->finalizeFluxData();
848 void writeWellspecReport(const SimulatorTimer& timer)
const
850 const auto changedWells = this->
schedule_
851 .changed_wells(timer.reportStepNum(), this->initialStep());
853 const auto changedWellLists = this->
schedule_
854 .changedWellLists(timer.reportStepNum(), this->initialStep());
856 if (changedWells.empty() && !changedWellLists) {
860 this->outputModule_->outputWellspecReport(changedWells,
862 timer.reportStepNum(),
863 timer.simulationTimeElapsed(),
864 timer.currentDateTime());
867 void writeWellflowReport(
const SimulatorTimer& timer,
869 const int wellsRequest)
const
871 this->outputModule_->outputTimeStamp(
"WELLS",
872 timer.simulationTimeElapsed(),
873 timer.reportStepNum(),
874 timer.currentDateTime());
876 const auto wantConnData = wellsRequest > 1;
878 this->outputModule_->outputProdLog(simStep, wantConnData);
879 this->outputModule_->outputInjLog(simStep, wantConnData);
880 this->outputModule_->outputCumLog(simStep, wantConnData);
881 this->outputModule_->outputMSWLog(simStep);
884 int initialStep()
const
886 const auto& initConfig = this->eclState().cfg().init();
888 return initConfig.restartRequested()
889 ? initConfig.getRestartStep()
893 Simulator& simulator_;
894 std::unique_ptr<OutputModule> outputModule_;
895 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:1102
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:969
const data::Solution & globalCellData() const
Definition: CollectDataOnIORank.hpp:93
Definition: EclGenericWriter.hpp:65
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:533
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:152
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 std::optional< Inplace > &initialInPlace, const InterRegFlowMap &interRegFlows, SummaryState &summaryState, UDQState &udqState)
Definition: EclGenericWriter_impl.hpp:626
const Schedule & schedule_
Definition: EclGenericWriter.hpp:155
SimulatorReport simulation_report_
Definition: EclGenericWriter.hpp:165
SimulatorReportSingle sub_step_report_
Definition: EclGenericWriter.hpp:164
const Dune::CartesianIndexMapper< GetPropType< TypeTag, Properties::Grid > > & cartMapper_
Definition: EclGenericWriter.hpp:161
std::unique_ptr< EclipseIO > eclIO_
Definition: EclGenericWriter.hpp:157
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:115
OutputModule & mutableOutputModule() const
Definition: EclWriter.hpp:716
const OutputModule & outputModule() const
Definition: EclWriter.hpp:713
void evalSummaryState(bool isSubStep)
collect and pass data and pass it to eclIO writer
Definition: EclWriter.hpp:203
static void registerParameters()
Definition: EclWriter.hpp:138
void serializeOp(Serializer &serializer)
Definition: EclWriter.hpp:723
void writeInitialFIPReport()
Writes the initial FIP report as configured in RPTSOL.
Definition: EclWriter.hpp:351
void beginRestart()
Definition: EclWriter.hpp:522
EclWriter(Simulator &simulator)
Definition: EclWriter.hpp:152
void writeOutput(data::Solution &&localCellData, bool isSubStep)
Definition: EclWriter.hpp:438
void writeReports(const SimulatorTimer &timer)
Definition: EclWriter.hpp:393
void endRestart()
Definition: EclWriter.hpp:699
Scalar restartTimeStepSize() const
Definition: EclWriter.hpp:719
~EclWriter()
Definition: EclWriter.hpp:192
const EquilGrid & globalGrid() const
Definition: EclWriter.hpp:195
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:43
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