EclWriter.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_ECL_WRITER_HPP
29#define OPM_ECL_WRITER_HPP
30
31#include <dune/grid/common/partitionset.hh>
32
33#include <opm/common/TimingMacros.hpp> // OPM_TIMEBLOCK
34#include <opm/common/OpmLog/OpmLog.hpp>
35#include <opm/input/eclipse/Schedule/RPTConfig.hpp>
36
37#include <opm/input/eclipse/Units/UnitSystem.hpp>
38#include <opm/input/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp>
39
40#include <opm/output/eclipse/Inplace.hpp>
41#include <opm/output/eclipse/RestartValue.hpp>
42
44#include <opm/models/blackoil/blackoilproperties.hh> // Properties::EnableMech, EnableSolvent
45#include <opm/models/common/multiphasebaseproperties.hh> // Properties::FluidSystem
46
55
56#include <boost/date_time/posix_time/posix_time.hpp>
57
58#include <limits>
59#include <map>
60#include <memory>
61#include <optional>
62#include <stdexcept>
63#include <string>
64#include <utility>
65#include <vector>
66
67namespace Opm::Parameters {
68
69// If available, write the ECL output in a non-blocking manner
70struct EnableAsyncEclOutput { static constexpr bool value = true; };
71
72// By default, use single precision for the ECL formated results
73struct EclOutputDoublePrecision { static constexpr bool value = false; };
74
75// Write all solutions for visualization, not just the ones for the
76// report steps...
77struct EnableWriteAllSolutions { static constexpr bool value = false; };
78
79// Write ESMRY file for fast loading of summary data
80struct EnableEsmry { static constexpr bool value = true; };
81
82} // namespace Opm::Parameters
83
84namespace Opm::Action {
85 class State;
86} // namespace Opm::Action
87
88namespace Opm {
89 class EclipseIO;
90 class UDQState;
91} // namespace Opm
92
93namespace Opm {
109template <class TypeTag, class OutputModule>
110class EclWriter : public EclGenericWriter<GetPropType<TypeTag, Properties::Grid>,
111 GetPropType<TypeTag, Properties::EquilGrid>,
112 GetPropType<TypeTag, Properties::GridView>,
113 GetPropType<TypeTag, Properties::ElementMapper>,
114 GetPropType<TypeTag, Properties::Scalar>>
115{
124 using Element = typename GridView::template Codim<0>::Entity;
126 using ElementIterator = typename GridView::template Codim<0>::Iterator;
128
129 typedef Dune::MultipleCodimMultipleGeomTypeMapper< GridView > VertexMapper;
130
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
136public:
137
138 static void registerParameters()
139 {
140 OutputModule::registerParameters();
141
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.");
147 }
148
149 // The Simulator object should preferably have been const - the
150 // only reason that is not the case is due to the SummaryState
151 // object owned deep down by the vanguard.
152 explicit EclWriter(Simulator& simulator)
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()
159 : nullptr),
160 simulator.vanguard().gridView(),
161 simulator.vanguard().cartesianIndexMapper(),
162 ((simulator.vanguard().grid().comm().rank() == 0)
163 ? &simulator.vanguard().equilCartesianIndexMapper()
164 : nullptr),
165 Parameters::Get<Parameters::EnableAsyncEclOutput>(),
166 Parameters::Get<Parameters::EnableEsmry>())
167 , simulator_(simulator)
168 {
169#if HAVE_MPI
170 if (this->simulator_.vanguard().grid().comm().size() > 1) {
171 auto smryCfg = (this->simulator_.vanguard().grid().comm().rank() == 0)
172 ? this->eclIO_->finalSummaryConfig()
173 : SummaryConfig{};
174
175 eclBroadcast(this->simulator_.vanguard().grid().comm(), smryCfg);
176
177 this->outputModule_ = std::make_unique<OutputModule>
178 (simulator, smryCfg, this->collectOnIORank_);
179 }
180 else
181#endif
182 {
183 this->outputModule_ = std::make_unique<OutputModule>
184 (simulator, this->eclIO_->finalSummaryConfig(), this->collectOnIORank_);
185 }
186
187 this->rank_ = this->simulator_.vanguard().grid().comm().rank();
188
189 this->simulator_.vanguard().eclState().computeFipRegionStatistics();
190 }
191
193 {}
194
195 const EquilGrid& globalGrid() const
196 {
197 return simulator_.vanguard().equilGrid();
198 }
199
203 void evalSummaryState(bool isSubStep)
204 {
205 OPM_TIMEBLOCK(evalSummaryState);
206 const int reportStepNum = simulator_.episodeIndex() + 1;
207
208 /*
209 The summary data is not evaluated for timestep 0, that is
210 implemented with a:
211
212 if (time_step == 0)
213 return;
214
215 check somewhere in the summary code. When the summary code was
216 split in separate methods Summary::eval() and
217 Summary::add_timestep() it was necessary to pull this test out
218 here to ensure that the well and group related keywords in the
219 restart file, like XWEL and XGRP were "correct" also in the
220 initial report step.
221
222 "Correct" in this context means unchanged behavior, might very
223 well be more correct to actually remove this if test.
224 */
225
226 if (reportStepNum == 0)
227 return;
228
229 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
230 const Scalar totalCpuTime =
231 simulator_.executionTimer().realTimeElapsed() +
232 simulator_.setupTimer().realTimeElapsed() +
233 simulator_.vanguard().setupTime();
234
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);
239
240 const auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
241 const auto localWellTestState = simulator_.problem().wellModel().wellTestState();
242 this->prepareLocalCellData(isSubStep, reportStepNum);
243
244 if (this->outputModule_->needInterfaceFluxes(isSubStep)) {
245 this->captureLocalFluxData();
246 }
247
248 if (this->collectOnIORank_.isParallel()) {
250
251 std::map<std::pair<std::string,int>,double> dummy;
252 this->collectOnIORank_.collect({},
253 outputModule_->getBlockData(),
254 dummy,
255 localWellData,
256 localWBP,
257 localGroupAndNetworkData,
258 localAquiferData,
259 localWellTestState,
260 this->outputModule_->getInterRegFlows(),
261 {},
262 {});
263
264 if (this->collectOnIORank_.isIORank()) {
265 auto& iregFlows = this->collectOnIORank_.globalInterRegFlows();
266
267 if (! iregFlows.readIsConsistent()) {
268 throw std::runtime_error {
269 "Inconsistent inter-region flow "
270 "region set names in parallel"
271 };
272 }
273
274 iregFlows.compress();
275 }
276
277 OPM_END_PARALLEL_TRY_CATCH("Collect to I/O rank: ",
278 this->simulator_.vanguard().grid().comm());
279 }
280
281
282 std::map<std::string, double> miscSummaryData;
283 std::map<std::string, std::vector<double>> regionData;
284 Inplace inplace;
285
286 {
287 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
288
289 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
290
291 if (this->collectOnIORank_.isIORank()){
292 inplace_ = inplace;
293 }
294 }
295
296 // Add TCPU
297 if (totalCpuTime != 0.0) {
298 miscSummaryData["TCPU"] = totalCpuTime;
299 }
301 miscSummaryData["NEWTON"] = this->sub_step_report_.total_newton_iterations;
302 }
304 miscSummaryData["MLINEARS"] = this->sub_step_report_.total_linear_iterations;
305 }
307 miscSummaryData["NLINEARS"] = static_cast<float>(this->sub_step_report_.total_linear_iterations) / this->sub_step_report_.total_newton_iterations;
308 }
309 if (this->sub_step_report_.min_linear_iterations != std::numeric_limits<unsigned int>::max()) {
310 miscSummaryData["NLINSMIN"] = this->sub_step_report_.min_linear_iterations;
311 }
313 miscSummaryData["NLINSMAX"] = this->sub_step_report_.max_linear_iterations;
314 }
316 miscSummaryData["MSUMLINS"] = this->simulation_report_.success.total_linear_iterations;
317 }
319 miscSummaryData["MSUMNEWT"] = this->simulation_report_.success.total_newton_iterations;
320 }
321
322 {
323 OPM_TIMEBLOCK(evalSummary);
324
325 const auto& blockData = this->collectOnIORank_.isParallel()
327 : this->outputModule_->getBlockData();
328
329 const auto& interRegFlows = this->collectOnIORank_.isParallel()
331 : this->outputModule_->getInterRegFlows();
332
333 this->evalSummary(reportStepNum,
334 curTime,
335 localWellData,
336 localWBP,
337 localGroupAndNetworkData,
338 localAquiferData,
339 blockData,
340 miscSummaryData,
341 regionData,
342 inplace,
343 this->outputModule_->initialInplace(),
344 interRegFlows,
345 this->summaryState(),
346 this->udqState());
347 }
348 }
349
352 {
353 const auto& gridView = simulator_.vanguard().gridView();
354 const int num_interior = detail::
356
357 this->outputModule_->
358 allocBuffers(num_interior, 0, false, false, /*isRestart*/ false);
359
360#ifdef _OPENMP
361#pragma omp parallel for
362#endif
363 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
364 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
365 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
366
367 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
368 }
369
370 // We always calculate the initial fip values as it may be used by various
371 // keywords in the Schedule, e.g. FIP=2 in RPTSCHED but no FIP in RPTSOL
372 outputModule_->calc_initial_inplace(simulator_.gridView().comm());
373
374 // check if RPTSOL entry has FIP output
375 const auto& fip = simulator_.vanguard().eclState().getEclipseConfig().fip();
376 if (fip.output(FIPConfig::OutputField::FIELD) ||
377 fip.output(FIPConfig::OutputField::RESV))
378 {
379 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
380 boost::posix_time::ptime start_time =
381 boost::posix_time::from_time_t(simulator_.vanguard().schedule().getStartTime());
382
383 if (this->collectOnIORank_.isIORank()) {
384 inplace_ = outputModule_->initialInplace().value();
385 outputModule_->outputFipAndResvLog(inplace_, 0, 0.0, start_time,
386 false, simulator_.gridView().comm());
387 }
388 }
389
390 outputModule_->outputFipAndResvLogToCSV(0, false, simulator_.gridView().comm());
391 }
392
393 void writeReports(const SimulatorTimer& timer)
394 {
395 if (! this->collectOnIORank_.isIORank()) {
396 return;
397 }
398
399 // SimulatorTimer::reportStepNum() is the simulator's zero-based
400 // "episode index". This is generally the index value needed to
401 // look up objects in the Schedule container. That said, function
402 // writeReports() is invoked at the *beginning* of a report
403 // step/episode which means we typically need the objects from the
404 // *previous* report step/episode. We therefore need special case
405 // handling for reportStepNum() == 0 in base runs and
406 // reportStepNum() <= restart step in restarted runs.
407 const auto firstStep = this->initialStep();
408 const auto simStep =
409 std::max(timer.reportStepNum() - 1, firstStep);
410
411 const auto& rpt = this->schedule_[simStep].rpt_config();
412
413 if (rpt.contains("WELSPECS") && (rpt.at("WELSPECS") > 0)) {
414 // Requesting a well specification report is valid at all times,
415 // including reportStepNum() == initialStep().
416 this->writeWellspecReport(timer);
417 }
418
419 if (timer.reportStepNum() == firstStep) {
420 // No dynamic flows at the beginning of the initialStep().
421 return;
422 }
423
424 if (rpt.contains("WELLS") && rpt.at("WELLS") > 0) {
425 this->writeWellflowReport(timer, simStep, rpt.at("WELLS"));
426 }
427
428 this->outputModule_->outputFipAndResvLog(this->inplace_,
429 timer.reportStepNum(),
430 timer.simulationTimeElapsed(),
431 timer.currentDateTime(),
432 /* isSubstep = */ false,
433 simulator_.gridView().comm());
434
435 OpmLog::note(""); // Blank line after all reports.
436 }
437
438 void writeOutput(data::Solution&& localCellData, bool isSubStep)
439 {
440 OPM_TIMEBLOCK(writeOutput);
441
442 const int reportStepNum = simulator_.episodeIndex() + 1;
443 this->prepareLocalCellData(isSubStep, reportStepNum);
444 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
445
446 // output using eclWriter if enabled
447 auto localWellData = simulator_.problem().wellModel().wellData();
448 auto localGroupAndNetworkData = simulator_.problem().wellModel()
449 .groupAndNetworkData(reportStepNum);
450
451 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
452 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
453
454 const bool isFlowsn = this->outputModule_->getFlows().hasFlowsn();
455 auto flowsn = this->outputModule_->getFlows().getFlowsn();
456
457 const bool isFloresn = this->outputModule_->getFlows().hasFloresn();
458 auto floresn = this->outputModule_->getFlows().getFloresn();
459
460 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
461
462 if (localCellData.empty()) {
463 this->outputModule_->assignToSolution(localCellData);
464 }
465
466 // Add cell data to perforations for RFT output
467 this->outputModule_->addRftDataToWells(localWellData,
468 reportStepNum,
469 simulator_.gridView().comm());
470 }
471
472 if (this->collectOnIORank_.isParallel() ||
473 this->collectOnIORank_.doesNeedReordering())
474 {
475 // Note: We don't need WBP (well-block averaged pressures) or
476 // inter-region flow rate values in order to create restart file
477 // output. There's consequently no need to collect those
478 // properties on the I/O rank.
479
480 this->collectOnIORank_.collect(localCellData,
481 this->outputModule_->getBlockData(),
482 this->outputModule_->getExtraBlockData(),
483 localWellData,
484 /* wbpData = */ {},
485 localGroupAndNetworkData,
486 localAquiferData,
487 localWellTestState,
488 /* interRegFlows = */ {},
489 flowsn,
490 floresn);
491 if (this->collectOnIORank_.isIORank()) {
492 this->outputModule_->assignGlobalFieldsToSolution(this->collectOnIORank_.globalCellData());
493 }
494 } else {
495 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
496 }
497
498 if (this->collectOnIORank_.isIORank()) {
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();
504 }
505 this->doWriteOutput(reportStepNum, timeStepIdx, isSubStep,
506 std::move(localCellData),
507 std::move(localWellData),
508 std::move(localGroupAndNetworkData),
509 std::move(localAquiferData),
510 std::move(localWellTestState),
511 this->actionState(),
512 this->udqState(),
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));
519 }
520 }
521
523 {
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");
531
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},
538
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()},
543
544 {"SGMAX", UnitSystem::measure::identity, enableNonWettingHysteresis && oilActive && gasActive},
545 {"SHMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && gasActive},
546
547 {"SOMAX", UnitSystem::measure::identity,
548 (enableNonWettingHysteresis && oilActive && waterActive)
549 || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
550
551 {"SOMIN", UnitSystem::measure::identity, enablePCHysteresis && oilActive && gasActive},
552 {"SWHY1", UnitSystem::measure::identity, enablePCHysteresis && oilActive && waterActive},
553 {"SWMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && waterActive},
554
555 {"PPCW", UnitSystem::measure::pressure, enableSwatinit},
556 };
557
558 {
559 const auto& tracers = simulator_.vanguard().eclState().tracer();
560
561 for (const auto& tracer : tracers) {
562 const auto enableSolTracer =
563 ((tracer.phase == Phase::GAS) && FluidSystem::enableDissolvedGas()) ||
564 ((tracer.phase == Phase::OIL) && FluidSystem::enableVaporizedOil());
565
566 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity, true);
567 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
568 }
569 }
570
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()},
575 };
576
577 const auto& gridView = this->simulator_.vanguard().gridView();
578 const auto numElements = gridView.size(/*codim=*/0);
579
580 // Try to load restart step 0 to calculate initial FIP
581 {
582 this->outputModule_->allocBuffers(numElements,
583 0,
584 /*isSubStep = */false,
585 /*log = */ false,
586 /*isRestart = */true);
587
588 const auto restartSolution =
590 solutionKeys, gridView.comm(), 0);
591
592 if (!restartSolution.empty()) {
593 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
594 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
595 this->outputModule_->setRestart(restartSolution, elemIdx, globalIdx);
596 }
597
598 this->simulator_.problem().readSolutionFromOutputModule(0, true);
599 this->simulator_.problem().temperatureModel().init();
600 ElementContext elemCtx(this->simulator_);
601 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
602 elemCtx.updatePrimaryStencil(elem);
603 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
604
605 this->outputModule_->updateFluidInPlace(elemCtx);
606 }
607
608 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
609 }
610 }
611
612 {
613 // The episodeIndex is rewound one step back before calling
614 // beginRestart() and cannot be used here. We just ask the
615 // initconfig directly to be sure that we use the correct index.
616 const auto restartStepIdx = this->simulator_.vanguard()
617 .eclState().getInitConfig().getRestartStep();
618
619 this->outputModule_->allocBuffers(numElements,
620 restartStepIdx,
621 /*isSubStep = */false,
622 /*log = */ false,
623 /*isRestart = */true);
624 }
625
626 {
627 const auto restartValues =
628 loadParallelRestart(this->eclIO_.get(),
629 this->actionState(),
630 this->summaryState(),
631 solutionKeys, extraKeys, gridView.comm());
632
633 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
634 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
635 this->outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
636 }
637
638 auto& tracer_model = simulator_.problem().tracerModel();
639 for (int tracer_index = 0; tracer_index < tracer_model.numTracers(); ++tracer_index) {
640 // Free tracers
641 {
642 const auto& free_tracer_name = tracer_model.fname(tracer_index);
643 const auto& free_tracer_solution = restartValues.solution
644 .template data<double>(free_tracer_name);
645
646 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
647 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
648 tracer_model.setFreeTracerConcentration
649 (tracer_index, elemIdx, free_tracer_solution[globalIdx]);
650 }
651 }
652
653 // Solution tracer (only if DISGAS/VAPOIL are active for gas/oil tracers)
654 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
655 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil()))
656 {
657 tracer_model.setEnableSolTracers(tracer_index, true);
658
659 const auto& sol_tracer_name = tracer_model.sname(tracer_index);
660 const auto& sol_tracer_solution = restartValues.solution
661 .template data<double>(sol_tracer_name);
662
663 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
664 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
665 tracer_model.setSolTracerConcentration
666 (tracer_index, elemIdx, sol_tracer_solution[globalIdx]);
667 }
668 }
669 else {
670 tracer_model.setEnableSolTracers(tracer_index, false);
671
672 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
673 tracer_model.setSolTracerConcentration(tracer_index, elemIdx, 0.0);
674 }
675 }
676 }
677
678 if (inputThpres.active()) {
679 const_cast<Simulator&>(this->simulator_)
680 .problem().thresholdPressure()
681 .setFromRestart(restartValues.getExtra("THRESHPR"));
682 }
683
684 restartTimeStepSize_ = restartValues.getExtra("OPMEXTRA")[0];
685 if (restartTimeStepSize_ <= 0) {
686 restartTimeStepSize_ = std::numeric_limits<double>::max();
687 }
688
689 // Initialize the well model from restart values
690 this->simulator_.problem().wellModel()
691 .initFromRestartFile(restartValues);
692
693 if (!restartValues.aquifer.empty()) {
694 this->simulator_.problem().mutableAquiferModel()
695 .initFromRestart(restartValues.aquifer);
696 }
697 }
698 }
699
701 {
702 // Calculate initial in-place volumes.
703 // Does nothing if they have already been calculated,
704 // e.g. from restart data at T=0.
705 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
706
707 if (this->collectOnIORank_.isIORank()) {
708 if (this->outputModule_->initialInplace().has_value()) {
709 this->inplace_ = this->outputModule_->initialInplace().value();
710 }
711 }
712 }
713
714 const OutputModule& outputModule() const
715 { return *outputModule_; }
716
717 OutputModule& mutableOutputModule() const
718 { return *outputModule_; }
719
720 Scalar restartTimeStepSize() const
721 { return restartTimeStepSize_; }
722
723 template <class Serializer>
724 void serializeOp(Serializer& serializer)
725 {
726 serializer(*outputModule_);
727 }
728
729private:
730 static bool enableEclOutput_()
731 {
732 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
733 return enable;
734 }
735
736 const EclipseState& eclState() const
737 { return simulator_.vanguard().eclState(); }
738
739 SummaryState& summaryState()
740 { return simulator_.vanguard().summaryState(); }
741
742 Action::State& actionState()
743 { return simulator_.vanguard().actionState(); }
744
745 UDQState& udqState()
746 { return simulator_.vanguard().udqState(); }
747
748 const Schedule& schedule() const
749 { return simulator_.vanguard().schedule(); }
750
751 void prepareLocalCellData(const bool isSubStep,
752 const int reportStepNum)
753 {
754 OPM_TIMEBLOCK(prepareLocalCellData);
755
756 if (this->outputModule_->localDataValid()) {
757 return;
758 }
759
760 const auto& gridView = simulator_.vanguard().gridView();
761 const bool log = this->collectOnIORank_.isIORank();
762
763 const int num_interior = detail::
765 this->outputModule_->
766 allocBuffers(num_interior, reportStepNum,
767 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
768 log, /*isRestart*/ false);
769
770 ElementContext elemCtx(simulator_);
771
773
774 {
775 OPM_TIMEBLOCK(prepareCellBasedData);
776
777 this->outputModule_->prepareDensityAccumulation();
778 this->outputModule_->setupExtractors(isSubStep, reportStepNum);
779 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
780 elemCtx.updatePrimaryStencil(elem);
781 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
782
783 this->outputModule_->processElement(elemCtx);
784 this->outputModule_->processElementBlockData(elemCtx);
785 }
786 this->outputModule_->clearExtractors();
787
788 this->outputModule_->accumulateDensityParallel();
789 }
790
791 {
792 OPM_TIMEBLOCK(prepareFluidInPlace);
793
794#ifdef _OPENMP
795#pragma omp parallel for
796#endif
797 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
798 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
799 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
800
801 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
802 }
803 }
804
805 this->outputModule_->validateLocalData();
806
807 OPM_END_PARALLEL_TRY_CATCH("EclWriter::prepareLocalCellData() failed: ",
808 this->simulator_.vanguard().grid().comm());
809 }
810
811 void captureLocalFluxData()
812 {
813 OPM_TIMEBLOCK(captureLocalData);
814
815 const auto& gridView = this->simulator_.vanguard().gridView();
816 const auto timeIdx = 0u;
817
818 auto elemCtx = ElementContext { this->simulator_ };
819
820 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
821 const auto activeIndex = [&elemMapper](const Element& e)
822 {
823 return elemMapper.index(e);
824 };
825
826 const auto cartesianIndex = [this](const int elemIndex)
827 {
828 return this->cartMapper_.cartesianIndex(elemIndex);
829 };
830
831 this->outputModule_->initializeFluxData();
832
834
835 for (const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
836 elemCtx.updateStencil(elem);
837 elemCtx.updateIntensiveQuantities(timeIdx);
838 elemCtx.updateExtensiveQuantities(timeIdx);
839
840 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
841 }
842
843 OPM_END_PARALLEL_TRY_CATCH("EclWriter::captureLocalFluxData() failed: ",
844 this->simulator_.vanguard().grid().comm())
845
846 this->outputModule_->finalizeFluxData();
847 }
848
849 void writeWellspecReport(const SimulatorTimer& timer) const
850 {
851 const auto changedWells = this->schedule_
852 .changed_wells(timer.reportStepNum(), this->initialStep());
853
854 const auto changedWellLists = this->schedule_
855 .changedWellLists(timer.reportStepNum(), this->initialStep());
856
857 if (changedWells.empty() && !changedWellLists) {
858 return;
859 }
860
861 this->outputModule_->outputWellspecReport(changedWells,
862 changedWellLists,
863 timer.reportStepNum(),
864 timer.simulationTimeElapsed(),
865 timer.currentDateTime());
866 }
867
868 void writeWellflowReport(const SimulatorTimer& timer,
869 const int simStep,
870 const int wellsRequest) const
871 {
872 this->outputModule_->outputTimeStamp("WELLS",
873 timer.simulationTimeElapsed(),
874 timer.reportStepNum(),
875 timer.currentDateTime());
876
877 const auto wantConnData = wellsRequest > 1;
878
879 this->outputModule_->outputProdLog(simStep, wantConnData);
880 this->outputModule_->outputInjLog(simStep, wantConnData);
881 this->outputModule_->outputCumLog(simStep, wantConnData);
882 this->outputModule_->outputMSWLog(simStep);
883 }
884
885 int initialStep() const
886 {
887 const auto& initConfig = this->eclState().cfg().init();
888
889 return initConfig.restartRequested()
890 ? initConfig.getRestartStep()
891 : 0;
892 }
893
894 Simulator& simulator_;
895 std::unique_ptr<OutputModule> outputModule_;
896 Scalar restartTimeStepSize_;
897 int rank_ ;
898 Inplace inplace_;
899};
900
901} // namespace Opm
902
903#endif // OPM_ECL_WRITER_HPP
#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:547
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 > > &regionData, const Inplace &inplace, const std::optional< Inplace > &initialInPlace, const InterRegFlowMap &interRegFlows, SummaryState &summaryState, UDQState &udqState)
Definition: EclGenericWriter_impl.hpp:655
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:115
OutputModule & mutableOutputModule() const
Definition: EclWriter.hpp:717
const OutputModule & outputModule() const
Definition: EclWriter.hpp:714
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:724
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:700
Scalar restartTimeStepSize() const
Definition: EclWriter.hpp:720
~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