28#ifndef OPM_ECL_WRITER_HPP
29#define OPM_ECL_WRITER_HPP
31#include <dune/grid/common/partitionset.hh>
33#include <opm/common/OpmLog/OpmLog.hpp>
34#include <opm/input/eclipse/Schedule/RPTConfig.hpp>
36#include <opm/input/eclipse/Units/UnitSystem.hpp>
38#include <opm/output/eclipse/RestartValue.hpp>
76namespace Action {
class State; }
95template <
class TypeTag>
97 GetPropType<TypeTag, Properties::EquilGrid>,
98 GetPropType<TypeTag, Properties::GridView>,
99 GetPropType<TypeTag, Properties::ElementMapper>,
100 GetPropType<TypeTag, Properties::Scalar>>
110 using Element =
typename GridView::template Codim<0>::Entity;
112 using ElementIterator =
typename GridView::template Codim<0>::Iterator;
115 typedef Dune::MultipleCodimMultipleGeomTypeMapper< GridView > VertexMapper;
117 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
118 enum { enableMech = getPropValue<TypeTag, Properties::EnableMech>() };
119 enum { enableTemperature = getPropValue<TypeTag, Properties::EnableTemperature>() };
120 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
128 Parameters::Register<Parameters::EnableAsyncEclOutput>
129 (
"Write the ECL-formated results in a non-blocking way "
130 "(i.e., using a separate thread).");
131 Parameters::Register<Parameters::EnableEsmry>
132 (
"Write ESMRY file for fast loading of summary data.");
139 :
BaseType(simulator.vanguard().schedule(),
140 simulator.vanguard().eclState(),
141 simulator.vanguard().summaryConfig(),
142 simulator.vanguard().grid(),
143 ((simulator.vanguard().grid().comm().rank() == 0)
144 ? &simulator.vanguard().equilGrid()
146 simulator.vanguard().gridView(),
147 simulator.vanguard().cartesianIndexMapper(),
148 ((simulator.vanguard().grid().comm().rank() == 0)
149 ? &simulator.vanguard().equilCartesianIndexMapper()
151 Parameters::
Get<Parameters::EnableAsyncEclOutput>(),
152 Parameters::
Get<Parameters::EnableEsmry>())
153 , simulator_(simulator)
156 if (this->simulator_.vanguard().grid().comm().size() > 1) {
157 auto smryCfg = (this->simulator_.vanguard().grid().comm().rank() == 0)
158 ? this->
eclIO_->finalSummaryConfig()
161 eclBroadcast(this->simulator_.vanguard().grid().comm(), smryCfg);
163 this->outputModule_ = std::make_unique<OutputBlackOilModule<TypeTag>>
169 this->outputModule_ = std::make_unique<OutputBlackOilModule<TypeTag>>
173 this->rank_ = this->simulator_.vanguard().grid().comm().rank();
175 this->simulator_.vanguard().eclState().computeFipRegionStatistics();
183 return simulator_.vanguard().equilGrid();
192 const int reportStepNum = simulator_.episodeIndex() + 1;
212 if (reportStepNum == 0)
215 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
216 const Scalar totalCpuTime =
217 simulator_.executionTimer().realTimeElapsed() +
218 simulator_.setupTimer().realTimeElapsed() +
219 simulator_.vanguard().setupTime();
221 const auto localWellData = simulator_.problem().wellModel().wellData();
222 const auto localWBP = simulator_.problem().wellModel().wellBlockAveragePressures();
223 const auto localGroupAndNetworkData = simulator_.problem().wellModel()
224 .groupAndNetworkData(reportStepNum);
226 const auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
227 const auto localWellTestState = simulator_.problem().wellModel().wellTestState();
228 this->prepareLocalCellData(isSubStep, reportStepNum);
230 if (this->outputModule_->needInterfaceFluxes(isSubStep)) {
231 this->captureLocalFluxData();
238 outputModule_->getBlockData(),
241 localGroupAndNetworkData,
244 this->outputModule_->getInterRegFlows(),
251 if (! iregFlows.readIsConsistent()) {
252 throw std::runtime_error {
253 "Inconsistent inter-region flow "
254 "region set names in parallel"
262 this->simulator_.vanguard().grid().comm());
266 std::map<std::string, double> miscSummaryData;
267 std::map<std::string, std::vector<double>> regionData;
271 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
273 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
281 if (totalCpuTime != 0.0) {
282 miscSummaryData[
"TCPU"] = totalCpuTime;
311 : this->outputModule_->getBlockData();
315 : this->outputModule_->getInterRegFlows();
321 localGroupAndNetworkData,
327 this->outputModule_->initialInplace(),
329 this->summaryState(),
337 const auto& fip = simulator_.vanguard().eclState().getEclipseConfig().fip();
338 if (!fip.output(FIPConfig::OutputField::FIELD) &&
339 !fip.output(FIPConfig::OutputField::RESV)) {
343 const auto& gridView = simulator_.vanguard().gridView();
347 this->outputModule_->
348 allocBuffers(num_interior, 0,
false,
false,
false);
351#pragma omp parallel for
353 for (
int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
354 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, 0);
355 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
357 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
360 std::map<std::string, double> miscSummaryData;
361 std::map<std::string, std::vector<double>> regionData;
364 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
366 boost::posix_time::ptime start_time = boost::posix_time::from_time_t(simulator_.vanguard().schedule().getStartTime());
368 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
373 outputModule_->outputFipAndResvLog(inplace_, 0, 0.0, start_time,
374 false, simulator_.gridView().comm());
384 const auto& rpt = this->
schedule_[rstep-1].rpt_config.get();
385 if (rpt.contains(
"WELLS") && rpt.at(
"WELLS") > 0) {
387 outputModule_->outputProdLog(rstep-1);
388 outputModule_->outputInjLog(rstep-1);
389 outputModule_->outputCumLog(rstep-1);
404 const int reportStepNum = simulator_.episodeIndex() + 1;
405 this->prepareLocalCellData(isSubStep, reportStepNum);
406 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
409 auto localWellData = simulator_.problem().wellModel().wellData();
410 auto localGroupAndNetworkData = simulator_.problem().wellModel()
411 .groupAndNetworkData(reportStepNum);
413 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
414 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
416 const bool isFlowsn = this->outputModule_->hasFlowsn();
417 auto flowsn = this->outputModule_->getFlowsn();
419 const bool isFloresn = this->outputModule_->hasFloresn();
420 auto floresn = this->outputModule_->getFloresn();
423 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
425 if (localCellData.empty()) {
426 this->outputModule_->assignToSolution(localCellData);
430 this->outputModule_->addRftDataToWells(localWellData, reportStepNum);
434 this->collectOnIORank_.doesNeedReordering())
442 this->outputModule_->getBlockData(),
445 localGroupAndNetworkData,
455 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
459 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
460 const Scalar nextStepSize = simulator_.problem().nextTimeStepSize();
461 std::optional<int> timeStepIdx;
462 if (Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
463 timeStepIdx = simulator_.timeStepIndex();
466 std::move(localCellData),
467 std::move(localWellData),
468 std::move(localGroupAndNetworkData),
469 std::move(localAquiferData),
470 std::move(localWellTestState),
473 this->summaryState(),
474 this->simulator_.problem().thresholdPressure().getRestartVector(),
475 curTime, nextStepSize,
476 Parameters::Get<Parameters::EclOutputDoublePrecision>(),
477 isFlowsn, std::move(flowsn),
478 isFloresn, std::move(floresn));
484 bool enablePCHysteresis = simulator_.problem().materialLawManager()->enablePCHysteresis();
485 bool enableNonWettingHysteresis = simulator_.problem().materialLawManager()->enableNonWettingHysteresis();
486 bool enableWettingHysteresis = simulator_.problem().materialLawManager()->enableWettingHysteresis();
487 bool oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
488 bool gasActive = FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
489 bool waterActive = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
490 bool enableSwatinit = simulator_.vanguard().eclState().fieldProps().has_double(
"SWATINIT");
491 bool opm_rst_file = Parameters::Get<Parameters::EnableOpmRstFile>();
492 bool read_temp = enableEnergy || (opm_rst_file && enableTemperature);
493 std::vector<RestartKey> solutionKeys{
494 {
"PRESSURE", UnitSystem::measure::pressure},
495 {
"SWAT", UnitSystem::measure::identity,
static_cast<bool>(waterActive)},
496 {
"SGAS", UnitSystem::measure::identity,
static_cast<bool>(gasActive)},
497 {
"TEMP" , UnitSystem::measure::temperature, read_temp},
498 {
"SSOLVENT" , UnitSystem::measure::identity, enableSolvent},
499 {
"RS", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGas()},
500 {
"RV", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedOil()},
501 {
"RVW", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedWater()},
502 {
"SGMAX", UnitSystem::measure::identity, (enableNonWettingHysteresis && oilActive && gasActive)},
503 {
"SHMAX", UnitSystem::measure::identity, (enableWettingHysteresis && oilActive && gasActive)},
504 {
"SOMAX", UnitSystem::measure::identity, (enableNonWettingHysteresis && oilActive && waterActive) || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
505 {
"SOMIN", UnitSystem::measure::identity, (enablePCHysteresis && oilActive && gasActive)},
506 {
"SWHY1", UnitSystem::measure::identity, (enablePCHysteresis && oilActive && waterActive)},
507 {
"SWMAX", UnitSystem::measure::identity, (enableWettingHysteresis && oilActive && waterActive)},
508 {
"PPCW", UnitSystem::measure::pressure, enableSwatinit}
511 const auto& inputThpres = eclState().getSimulationConfig().getThresholdPressure();
512 std::vector<RestartKey> extraKeys = {{
"OPMEXTRA", UnitSystem::measure::identity,
false},
513 {
"THRESHPR", UnitSystem::measure::pressure, inputThpres.active()}};
516 const auto& tracers = simulator_.vanguard().eclState().tracer();
517 for (
const auto& tracer : tracers) {
518 bool enableSolTracer = (tracer.phase == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
519 (tracer.phase == Phase::OIL && FluidSystem::enableVaporizedOil());
520 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity,
true);
521 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
529 const auto& initconfig = simulator_.vanguard().eclState().getInitConfig();
530 int restartStepIdx = initconfig.getRestartStep();
532 const auto& gridView = simulator_.vanguard().gridView();
533 unsigned numElements = gridView.size(0);
534 outputModule_->allocBuffers(numElements, restartStepIdx,
false,
false,
true);
537 SummaryState& summaryState = simulator_.vanguard().summaryState();
538 Action::State& actionState = simulator_.vanguard().actionState();
540 gridView.grid().comm());
541 for (
unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
543 outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
546 auto& tracer_model = simulator_.problem().tracerModel();
547 for (
int tracer_index = 0; tracer_index < tracer_model.numTracers(); tracer_index++) {
549 const auto& free_tracer_name = tracer_model.fname(tracer_index);
550 const auto& free_tracer_solution = restartValues.solution.template data<double>(free_tracer_name);
551 for (
unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
553 tracer_model.setFreeTracerConcentration(tracer_index, globalIdx, free_tracer_solution[globalIdx]);
556 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
557 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil())) {
558 tracer_model.setEnableSolTracers(tracer_index,
true);
559 const auto& sol_tracer_name = tracer_model.sname(tracer_index);
560 const auto& sol_tracer_solution = restartValues.solution.template data<double>(sol_tracer_name);
561 for (
unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
563 tracer_model.setSolTracerConcentration(tracer_index, globalIdx, sol_tracer_solution[globalIdx]);
567 tracer_model.setEnableSolTracers(tracer_index,
false);
568 for (
unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
570 tracer_model.setSolTracerConcentration(tracer_index, globalIdx, 0.0);
575 if (inputThpres.active()) {
576 Simulator& mutableSimulator =
const_cast<Simulator&
>(simulator_);
577 auto& thpres = mutableSimulator.problem().thresholdPressure();
578 const auto& thpresValues = restartValues.getExtra(
"THRESHPR");
579 thpres.setFromRestart(thpresValues);
581 restartTimeStepSize_ = restartValues.getExtra(
"OPMEXTRA")[0];
582 if (restartTimeStepSize_ <= 0)
583 restartTimeStepSize_ = std::numeric_limits<double>::max();
586 simulator_.problem().wellModel().initFromRestartFile(restartValues);
588 if (!restartValues.aquifer.empty())
589 simulator_.problem().mutableAquiferModel().initFromRestart(restartValues.aquifer);
597 {
return *outputModule_; }
600 {
return *outputModule_; }
603 {
return restartTimeStepSize_; }
605 template <
class Serializer>
608 serializer(*outputModule_);
612 static bool enableEclOutput_()
614 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
618 const EclipseState& eclState()
const
619 {
return simulator_.vanguard().eclState(); }
621 SummaryState& summaryState()
622 {
return simulator_.vanguard().summaryState(); }
624 Action::State& actionState()
625 {
return simulator_.vanguard().actionState(); }
628 {
return simulator_.vanguard().udqState(); }
630 const Schedule& schedule()
const
631 {
return simulator_.vanguard().schedule(); }
633 void prepareLocalCellData(
const bool isSubStep,
634 const int reportStepNum)
636 OPM_TIMEBLOCK(prepareLocalCellData);
638 if (this->outputModule_->localDataValid()) {
642 const auto& gridView = simulator_.vanguard().gridView();
647 this->outputModule_->
648 allocBuffers(num_interior, reportStepNum,
649 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
652 ElementContext elemCtx(simulator_);
657 OPM_TIMEBLOCK(prepareCellBasedData);
659 this->outputModule_->prepareDensityAccumulation();
661 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
662 elemCtx.updatePrimaryStencil(elem);
663 elemCtx.updatePrimaryIntensiveQuantities(0);
665 this->outputModule_->processElement(elemCtx);
668 this->outputModule_->accumulateDensityParallel();
671 if constexpr (enableMech) {
672 if (simulator_.vanguard().eclState().runspec().mech()) {
673 OPM_TIMEBLOCK(prepareMechData);
674 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
675 elemCtx.updatePrimaryStencil(elem);
676 elemCtx.updatePrimaryIntensiveQuantities(0);
677 outputModule_->processElementMech(elemCtx);
682 if (! this->simulator_.model().linearizer().getFlowsInfo().empty()) {
683 OPM_TIMEBLOCK(prepareFlowsData);
684 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
685 elemCtx.updatePrimaryStencil(elem);
686 elemCtx.updatePrimaryIntensiveQuantities(0);
688 this->outputModule_->processElementFlows(elemCtx);
693 OPM_TIMEBLOCK(prepareBlockData);
694 for (
const auto& elem : elements(gridView, Dune::Partitions::interior)) {
695 elemCtx.updatePrimaryStencil(elem);
696 elemCtx.updatePrimaryIntensiveQuantities(0);
698 this->outputModule_->processElementBlockData(elemCtx);
703 OPM_TIMEBLOCK(prepareFluidInPlace);
706#pragma omp parallel for
708 for (
int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
709 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, 0);
710 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
712 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
716 this->outputModule_->validateLocalData();
719 this->simulator_.vanguard().grid().comm());
722 void captureLocalFluxData()
724 OPM_TIMEBLOCK(captureLocalData);
726 const auto& gridView = this->simulator_.vanguard().gridView();
727 const auto timeIdx = 0u;
729 auto elemCtx = ElementContext { this->simulator_ };
731 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
732 const auto activeIndex = [&elemMapper](
const Element& e)
734 return elemMapper.index(e);
737 const auto cartesianIndex = [
this](
const int elemIndex)
739 return this->
cartMapper_.cartesianIndex(elemIndex);
742 this->outputModule_->initializeFluxData();
746 for (
const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
747 elemCtx.updateStencil(elem);
748 elemCtx.updateIntensiveQuantities(timeIdx);
749 elemCtx.updateExtensiveQuantities(timeIdx);
751 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
755 this->simulator_.vanguard().grid().comm())
757 this->outputModule_->finalizeFluxData();
760 Simulator& simulator_;
761 std::unique_ptr<OutputBlackOilModule<TypeTag> > outputModule_;
762 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
void collect(const data::Solution &localCellData, const std::map< std::pair< std::string, int >, double > &localBlockData, 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
int localIdxToGlobalIdx(unsigned localIdx) const
Definition: CollectDataOnIORank_impl.hpp:1091
InterRegFlowMap & globalInterRegFlows()
Definition: CollectDataOnIORank.hpp:113
bool isParallel() const
Definition: CollectDataOnIORank.hpp:128
bool isIORank() const
Definition: CollectDataOnIORank.hpp:125
const std::map< std::pair< std::string, int >, double > & globalBlockData() const
Definition: CollectDataOnIORank.hpp:89
const data::Solution & globalCellData() const
Definition: CollectDataOnIORank.hpp:92
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:531
CollectDataOnIORankType collectOnIORank_
Definition: EclGenericWriter.hpp:152
const Schedule & schedule_
Definition: EclGenericWriter.hpp:155
SimulatorReport simulation_report_
Definition: EclGenericWriter.hpp:165
SimulatorReportSingle sub_step_report_
Definition: EclGenericWriter.hpp:164
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:624
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:101
void writeInitialFIPReport()
Writes the initial FIP report as configured in RPTSOL.
Definition: EclWriter.hpp:335
void endRestart()
Definition: EclWriter.hpp:593
const OutputBlackOilModule< TypeTag > & outputModule() const
Definition: EclWriter.hpp:596
Scalar restartTimeStepSize() const
Definition: EclWriter.hpp:602
~EclWriter()
Definition: EclWriter.hpp:178
void writeOutput(data::Solution &&localCellData, bool isSubStep)
Definition: EclWriter.hpp:400
EclWriter(Simulator &simulator)
Definition: EclWriter.hpp:138
static void registerParameters()
Definition: EclWriter.hpp:124
void writeReports(const SimulatorTimer &timer)
Definition: EclWriter.hpp:379
const EquilGrid & globalGrid() const
Definition: EclWriter.hpp:181
void beginRestart()
Definition: EclWriter.hpp:482
OutputBlackOilModule< TypeTag > & mutableOutputModule() const
Definition: EclWriter.hpp:599
void evalSummaryState(bool isSubStep)
collect and pass data and pass it to eclIO writer
Definition: EclWriter.hpp:189
void serializeOp(Serializer &serializer)
Definition: EclWriter.hpp:606
Output module for the results black oil model writing in ECL binary format.
Definition: OutputBlackoilModule.hpp:90
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: OutputBlackoilModule.hpp:182
virtual int reportStepNum() const
Current report step number. This might differ from currentStepNum in case of sub stepping.
Definition: SimulatorTimerInterface.hpp:50
Definition: SimulatorTimer.hpp:39
virtual boost::posix_time::ptime currentDateTime() const
Return the current time as a posix time object.
double simulationTimeElapsed() const override
Definition: blackoilnewtonmethodparams.hpp:31
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:185
std::size_t countLocalInteriorCellsGridView(const GridView &gridView)
Get the number of local interior cells in a grid view.
Definition: countGlobalCells.hpp:47
Definition: blackoilboundaryratevector.hh:37
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:235
Definition: EclWriter.hpp:63
static constexpr bool value
Definition: EclWriter.hpp:63
Definition: EclWriter.hpp:60
static constexpr bool value
Definition: EclWriter.hpp:60
Definition: EclWriter.hpp:57
static constexpr bool value
Definition: EclWriter.hpp:57
Definition: EclWriter.hpp:70
static constexpr bool value
Definition: EclWriter.hpp:70
Definition: EclWriter.hpp:67
static constexpr bool value
Definition: EclWriter.hpp:67
SimulatorReportSingle success
Definition: SimulatorReport.hpp:101
unsigned int min_linear_iterations
Definition: SimulatorReport.hpp:51
unsigned int total_newton_iterations
Definition: SimulatorReport.hpp:49
unsigned int max_linear_iterations
Definition: SimulatorReport.hpp:52
unsigned int total_linear_iterations
Definition: SimulatorReport.hpp:50