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 enum { enableGeochemistry = getPropValue<TypeTag, Properties::EnableGeochemistry>() };
136
137public:
138
139 static void registerParameters()
140 {
141 OutputModule::registerParameters();
142
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.");
148 }
149
150 // The Simulator object should preferably have been const - the
151 // only reason that is not the case is due to the SummaryState
152 // object owned deep down by the vanguard.
153 explicit EclWriter(Simulator& simulator)
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()
160 : nullptr),
161 simulator.vanguard().gridView(),
162 simulator.vanguard().cartesianIndexMapper(),
163 ((simulator.vanguard().grid().comm().rank() == 0)
164 ? &simulator.vanguard().equilCartesianIndexMapper()
165 : nullptr),
166 Parameters::Get<Parameters::EnableAsyncEclOutput>(),
167 Parameters::Get<Parameters::EnableEsmry>())
168 , simulator_(simulator)
169 {
170#if HAVE_MPI
171 if (this->simulator_.vanguard().grid().comm().size() > 1) {
172 auto smryCfg = (this->simulator_.vanguard().grid().comm().rank() == 0)
173 ? this->eclIO_->finalSummaryConfig()
174 : SummaryConfig{};
175
176 eclBroadcast(this->simulator_.vanguard().grid().comm(), smryCfg);
177
178 this->outputModule_ = std::make_unique<OutputModule>
179 (simulator, smryCfg, this->collectOnIORank_);
180 }
181 else
182#endif
183 {
184 this->outputModule_ = std::make_unique<OutputModule>
185 (simulator, this->eclIO_->finalSummaryConfig(), this->collectOnIORank_);
186 }
187
188 this->rank_ = this->simulator_.vanguard().grid().comm().rank();
189
190 this->simulator_.vanguard().eclState().computeFipRegionStatistics();
191 }
192
194 {}
195
196 const EquilGrid& globalGrid() const
197 {
198 return simulator_.vanguard().equilGrid();
199 }
200
204 void evalSummaryState(bool isSubStep)
205 {
206 OPM_TIMEBLOCK(evalSummaryState);
207 const int reportStepNum = simulator_.episodeIndex() + 1;
208
209 /*
210 The summary data is not evaluated for timestep 0, that is
211 implemented with a:
212
213 if (time_step == 0)
214 return;
215
216 check somewhere in the summary code. When the summary code was
217 split in separate methods Summary::eval() and
218 Summary::add_timestep() it was necessary to pull this test out
219 here to ensure that the well and group related keywords in the
220 restart file, like XWEL and XGRP were "correct" also in the
221 initial report step.
222
223 "Correct" in this context means unchanged behavior, might very
224 well be more correct to actually remove this if test.
225 */
226
227 if (reportStepNum == 0)
228 return;
229
230 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
231 const Scalar totalCpuTime =
232 simulator_.executionTimer().realTimeElapsed() +
233 simulator_.setupTimer().realTimeElapsed() +
234 simulator_.vanguard().setupTime();
235
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);
240
241 const auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
242 const auto localWellTestState = simulator_.problem().wellModel().wellTestState();
243 this->prepareLocalCellData(isSubStep, reportStepNum);
244
245 if (this->outputModule_->needInterfaceFluxes(isSubStep)) {
246 this->captureLocalFluxData();
247 }
248
249 if (this->collectOnIORank_.isParallel()) {
251
252 std::map<std::pair<std::string,int>,double> dummy;
253 this->collectOnIORank_.collect({},
254 outputModule_->getBlockData(),
255 dummy,
256 localWellData,
257 localWBP,
258 localGroupAndNetworkData,
259 localAquiferData,
260 localWellTestState,
261 this->outputModule_->getInterRegFlows(),
262 {},
263 {});
264
265 if (this->collectOnIORank_.isIORank()) {
266 auto& iregFlows = this->collectOnIORank_.globalInterRegFlows();
267
268 if (! iregFlows.readIsConsistent()) {
269 throw std::runtime_error {
270 "Inconsistent inter-region flow "
271 "region set names in parallel"
272 };
273 }
274
275 iregFlows.compress();
276 }
277
278 OPM_END_PARALLEL_TRY_CATCH("Collect to I/O rank: ",
279 this->simulator_.vanguard().grid().comm());
280 }
281
282
283 std::map<std::string, double> miscSummaryData;
284 std::map<std::string, std::vector<double>> regionData;
285 Inplace inplace;
286
287 {
288 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
289
290 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
291
292 if (this->collectOnIORank_.isIORank()){
293 inplace_ = inplace;
294 }
295 }
296
297 // Add TCPU
298 if (totalCpuTime != 0.0) {
299 miscSummaryData["TCPU"] = totalCpuTime;
300 }
302 miscSummaryData["NEWTON"] = this->sub_step_report_.total_newton_iterations;
303 }
305 miscSummaryData["MLINEARS"] = this->sub_step_report_.total_linear_iterations;
306 }
308 miscSummaryData["NLINEARS"] = static_cast<float>(this->sub_step_report_.total_linear_iterations) / this->sub_step_report_.total_newton_iterations;
309 }
310 if (this->sub_step_report_.min_linear_iterations != std::numeric_limits<unsigned int>::max()) {
311 miscSummaryData["NLINSMIN"] = this->sub_step_report_.min_linear_iterations;
312 }
314 miscSummaryData["NLINSMAX"] = this->sub_step_report_.max_linear_iterations;
315 }
317 miscSummaryData["MSUMLINS"] = this->simulation_report_.success.total_linear_iterations;
318 }
320 miscSummaryData["MSUMNEWT"] = this->simulation_report_.success.total_newton_iterations;
321 }
322
323 {
324 OPM_TIMEBLOCK(evalSummary);
325
326 const auto& blockData = this->collectOnIORank_.isParallel()
328 : this->outputModule_->getBlockData();
329
330 const auto& interRegFlows = this->collectOnIORank_.isParallel()
332 : this->outputModule_->getInterRegFlows();
333
334 this->evalSummary(reportStepNum,
335 curTime,
336 localWellData,
337 localWBP,
338 localGroupAndNetworkData,
339 localAquiferData,
340 blockData,
341 miscSummaryData,
342 regionData,
343 inplace,
344 this->outputModule_->initialInplace(),
345 interRegFlows,
346 this->summaryState(),
347 this->udqState());
348 }
349 }
350
353 {
354 const auto& gridView = simulator_.vanguard().gridView();
355 const int num_interior = detail::
357
358 this->outputModule_->
359 allocBuffers(num_interior, 0, false, false, /*isRestart*/ false);
360
361#ifdef _OPENMP
362#pragma omp parallel for
363#endif
364 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
365 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
366 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
367
368 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
369 }
370
371 // We always calculate the initial fip values as it may be used by various
372 // keywords in the Schedule, e.g. FIP=2 in RPTSCHED but no FIP in RPTSOL
373 outputModule_->calc_initial_inplace(simulator_.gridView().comm());
374
375 // check if RPTSOL entry has FIP output
376 const auto& fip = simulator_.vanguard().eclState().getEclipseConfig().fip();
377 if (fip.output(FIPConfig::OutputField::FIELD) ||
378 fip.output(FIPConfig::OutputField::RESV))
379 {
380 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
381
382 const auto start_time = boost::posix_time::
383 from_time_t(simulator_.vanguard().schedule().getStartTime());
384
385 if (this->collectOnIORank_.isIORank()) {
386 this->inplace_ = *this->outputModule_->initialInplace();
387
388 this->outputModule_->
389 outputFipAndResvLog(this->inplace_, 0, 0.0, start_time,
390 false, simulator_.gridView().comm());
391 }
392 }
393
394 outputModule_->outputFipAndResvLogToCSV(0, false, simulator_.gridView().comm());
395 }
396
397 void writeReports(const SimulatorTimer& timer)
398 {
399 if (! this->collectOnIORank_.isIORank()) {
400 return;
401 }
402
403 // SimulatorTimer::reportStepNum() is the simulator's zero-based
404 // "episode index". This is generally the index value needed to
405 // look up objects in the Schedule container. That said, function
406 // writeReports() is invoked at the *beginning* of a report
407 // step/episode which means we typically need the objects from the
408 // *previous* report step/episode. We therefore need special case
409 // handling for reportStepNum() == 0 in base runs and
410 // reportStepNum() <= restart step in restarted runs.
411 const auto firstStep = this->initialStep();
412 const auto simStep =
413 std::max(timer.reportStepNum() - 1, firstStep);
414
415 const auto& rpt = this->schedule_[simStep].rpt_config();
416
417 if (rpt.contains("WELSPECS") && (rpt.at("WELSPECS") > 0)) {
418 // Requesting a well specification report is valid at all times,
419 // including reportStepNum() == initialStep().
420 this->writeWellspecReport(timer);
421 }
422
423 if (timer.reportStepNum() == firstStep) {
424 // No dynamic flows at the beginning of the initialStep().
425 return;
426 }
427
428 if (rpt.contains("WELLS") && rpt.at("WELLS") > 0) {
429 this->writeWellflowReport(timer, simStep, rpt.at("WELLS"));
430 }
431
432 this->outputModule_->outputFipAndResvLog(this->inplace_,
433 timer.reportStepNum(),
434 timer.simulationTimeElapsed(),
435 timer.currentDateTime(),
436 /* isSubstep = */ false,
437 simulator_.gridView().comm());
438
439 OpmLog::note(""); // Blank line after all reports.
440 }
441
442 void writeOutput(data::Solution&& localCellData, bool isSubStep)
443 {
444 OPM_TIMEBLOCK(writeOutput);
445
446 const int reportStepNum = simulator_.episodeIndex() + 1;
447 this->prepareLocalCellData(isSubStep, reportStepNum);
448 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
449
450 // output using eclWriter if enabled
451 auto localWellData = simulator_.problem().wellModel().wellData();
452 auto localGroupAndNetworkData = simulator_.problem().wellModel()
453 .groupAndNetworkData(reportStepNum);
454
455 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
456 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
457
458 const bool isFlowsn = this->outputModule_->getFlows().hasFlowsn();
459 auto flowsn = this->outputModule_->getFlows().getFlowsn();
460
461 const bool isFloresn = this->outputModule_->getFlows().hasFloresn();
462 auto floresn = this->outputModule_->getFlows().getFloresn();
463
464 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
465
466 if (localCellData.empty()) {
467 this->outputModule_->assignToSolution(localCellData);
468 }
469
470 // Add cell data to perforations for RFT output
471 this->outputModule_->addRftDataToWells(localWellData,
472 reportStepNum,
473 simulator_.gridView().comm());
474 }
475
476 if (this->collectOnIORank_.isParallel() ||
477 this->collectOnIORank_.doesNeedReordering())
478 {
479 // Note: We don't need WBP (well-block averaged pressures) or
480 // inter-region flow rate values in order to create restart file
481 // output. There's consequently no need to collect those
482 // properties on the I/O rank.
483
484 this->collectOnIORank_.collect(localCellData,
485 this->outputModule_->getBlockData(),
486 this->outputModule_->getExtraBlockData(),
487 localWellData,
488 /* wbpData = */ {},
489 localGroupAndNetworkData,
490 localAquiferData,
491 localWellTestState,
492 /* interRegFlows = */ {},
493 flowsn,
494 floresn);
495 if (this->collectOnIORank_.isIORank()) {
496 this->outputModule_->assignGlobalFieldsToSolution(this->collectOnIORank_.globalCellData());
497 }
498 } else {
499 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
500 }
501
502 if (this->collectOnIORank_.isIORank()) {
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();
508 }
509 this->doWriteOutput(reportStepNum, timeStepIdx, isSubStep,
510 std::move(localCellData),
511 std::move(localWellData),
512 std::move(localGroupAndNetworkData),
513 std::move(localAquiferData),
514 std::move(localWellTestState),
515 this->actionState(),
516 this->udqState(),
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));
523 }
524 }
525
527 {
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");
535
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},
542
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()},
547
548 {"SGMAX", UnitSystem::measure::identity, enableNonWettingHysteresis && oilActive && gasActive},
549 {"SHMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && gasActive},
550
551 {"SOMAX", UnitSystem::measure::identity,
552 (enableNonWettingHysteresis && oilActive && waterActive)
553 || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
554
555 {"SOMIN", UnitSystem::measure::identity, enablePCHysteresis && oilActive && gasActive},
556 {"SWHY1", UnitSystem::measure::identity, enablePCHysteresis && oilActive && waterActive},
557 {"SWMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && waterActive},
558
559 {"PPCW", UnitSystem::measure::pressure, enableSwatinit},
560 };
561
562 {
563 const auto& tracers = simulator_.vanguard().eclState().tracer();
564
565 for (const auto& tracer : tracers) {
566 const auto enableSolTracer =
567 ((tracer.phase == Phase::GAS) && FluidSystem::enableDissolvedGas()) ||
568 ((tracer.phase == Phase::OIL) && FluidSystem::enableVaporizedOil());
569
570 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity, true);
571 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
572 }
573 }
574
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()},
579 };
580
581 const auto& gridView = this->simulator_.vanguard().gridView();
582 const auto numElements = gridView.size(/*codim=*/0);
583
584 // Try to load restart step 0 to calculate initial FIP
585 {
586 this->outputModule_->allocBuffers(numElements,
587 0,
588 /*isSubStep = */false,
589 /*log = */ false,
590 /*isRestart = */true);
591
592 const auto restartSolution =
594 solutionKeys, gridView.comm(), 0);
595
596 if (!restartSolution.empty()) {
597 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
598 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
599 this->outputModule_->setRestart(restartSolution, elemIdx, globalIdx);
600 }
601
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(/*timeIdx=*/0);
608
609 this->outputModule_->updateFluidInPlace(elemCtx);
610 }
611
612 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
613 }
614 }
615
616 {
617 // The episodeIndex is rewound one step back before calling
618 // beginRestart() and cannot be used here. We just ask the
619 // initconfig directly to be sure that we use the correct index.
620 const auto restartStepIdx = this->simulator_.vanguard()
621 .eclState().getInitConfig().getRestartStep();
622
623 this->outputModule_->allocBuffers(numElements,
624 restartStepIdx,
625 /*isSubStep = */false,
626 /*log = */ false,
627 /*isRestart = */true);
628 }
629
630 {
631 const auto restartValues =
632 loadParallelRestart(this->eclIO_.get(),
633 this->actionState(),
634 this->summaryState(),
635 solutionKeys, extraKeys, gridView.comm());
636
637 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
638 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
639 this->outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
640 }
641
642 auto& tracer_model = simulator_.problem().tracerModel();
643 for (int tracer_index = 0; tracer_index < tracer_model.numTracers(); ++tracer_index) {
644 // Free tracers
645 {
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);
649
650 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
651 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
652 tracer_model.setFreeTracerConcentration
653 (tracer_index, elemIdx, free_tracer_solution[globalIdx]);
654 }
655 }
656
657 // Solution tracer (only if DISGAS/VAPOIL are active for gas/oil tracers)
658 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
659 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil()))
660 {
661 tracer_model.setEnableSolTracers(tracer_index, true);
662
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);
666
667 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
668 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
669 tracer_model.setSolTracerConcentration
670 (tracer_index, elemIdx, sol_tracer_solution[globalIdx]);
671 }
672 }
673 else {
674 tracer_model.setEnableSolTracers(tracer_index, false);
675
676 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
677 tracer_model.setSolTracerConcentration(tracer_index, elemIdx, 0.0);
678 }
679 }
680 }
681
682 if (inputThpres.active()) {
683 const_cast<Simulator&>(this->simulator_)
684 .problem().thresholdPressure()
685 .setFromRestart(restartValues.getExtra("THRESHPR"));
686 }
687
688 restartTimeStepSize_ = restartValues.getExtra("OPMEXTRA")[0];
689 if (restartTimeStepSize_ <= 0) {
690 restartTimeStepSize_ = std::numeric_limits<double>::max();
691 }
692
693 // Initialize the well model from restart values
694 this->simulator_.problem().wellModel()
695 .initFromRestartFile(restartValues);
696
697 if (!restartValues.aquifer.empty()) {
698 this->simulator_.problem().mutableAquiferModel()
699 .initFromRestart(restartValues.aquifer);
700 }
701 }
702 }
703
705 {
706 // Calculate initial in-place volumes.
707 // Does nothing if they have already been calculated,
708 // e.g. from restart data at T=0.
709 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
710
711 if (this->collectOnIORank_.isIORank()) {
712 if (const auto* iip = this->outputModule_->initialInplace(); iip != nullptr) {
713 this->inplace_ = *iip;
714 }
715 }
716 }
717
718 const OutputModule& outputModule() const
719 { return *outputModule_; }
720
721 OutputModule& mutableOutputModule() const
722 { return *outputModule_; }
723
724 Scalar restartTimeStepSize() const
725 { return restartTimeStepSize_; }
726
727 template <class Serializer>
728 void serializeOp(Serializer& serializer)
729 {
730 serializer(*outputModule_);
731 }
732
733private:
734 static bool enableEclOutput_()
735 {
736 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
737 return enable;
738 }
739
740 const EclipseState& eclState() const
741 { return simulator_.vanguard().eclState(); }
742
743 SummaryState& summaryState()
744 { return simulator_.vanguard().summaryState(); }
745
746 Action::State& actionState()
747 { return simulator_.vanguard().actionState(); }
748
749 UDQState& udqState()
750 { return simulator_.vanguard().udqState(); }
751
752 const Schedule& schedule() const
753 { return simulator_.vanguard().schedule(); }
754
755 void prepareLocalCellData(const bool isSubStep,
756 const int reportStepNum)
757 {
758 OPM_TIMEBLOCK(prepareLocalCellData);
759
760 if (this->outputModule_->localDataValid()) {
761 return;
762 }
763
764 const auto& gridView = simulator_.vanguard().gridView();
765 const bool log = this->collectOnIORank_.isIORank();
766
767 const int num_interior = detail::
769 this->outputModule_->
770 allocBuffers(num_interior, reportStepNum,
771 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
772 log, /*isRestart*/ false);
773
774 ElementContext elemCtx(simulator_);
775
777
778 {
779 OPM_TIMEBLOCK(prepareCellBasedData);
780
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(/*timeIdx=*/0);
786
787 this->outputModule_->processElement(elemCtx);
788 this->outputModule_->processElementBlockData(elemCtx);
789 }
790 this->outputModule_->clearExtractors();
791
792 this->outputModule_->accumulateDensityParallel();
793 }
794
795 {
796 OPM_TIMEBLOCK(prepareFluidInPlace);
797
798#ifdef _OPENMP
799#pragma omp parallel for
800#endif
801 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
802 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
803 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
804
805 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
806 }
807 }
808
809 this->outputModule_->validateLocalData();
810
811 OPM_END_PARALLEL_TRY_CATCH("EclWriter::prepareLocalCellData() failed: ",
812 this->simulator_.vanguard().grid().comm());
813 }
814
815 void captureLocalFluxData()
816 {
817 OPM_TIMEBLOCK(captureLocalData);
818
819 const auto& gridView = this->simulator_.vanguard().gridView();
820 const auto timeIdx = 0u;
821
822 auto elemCtx = ElementContext { this->simulator_ };
823
824 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
825 const auto activeIndex = [&elemMapper](const Element& e)
826 {
827 return elemMapper.index(e);
828 };
829
830 const auto cartesianIndex = [this](const int elemIndex)
831 {
832 return this->cartMapper_.cartesianIndex(elemIndex);
833 };
834
835 this->outputModule_->initializeFluxData();
836
838
839 for (const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
840 elemCtx.updateStencil(elem);
841 elemCtx.updateIntensiveQuantities(timeIdx);
842 elemCtx.updateExtensiveQuantities(timeIdx);
843
844 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
845 }
846
847 OPM_END_PARALLEL_TRY_CATCH("EclWriter::captureLocalFluxData() failed: ",
848 this->simulator_.vanguard().grid().comm())
849
850 this->outputModule_->finalizeFluxData();
851 }
852
853 void writeWellspecReport(const SimulatorTimer& timer) const
854 {
855 const auto changedWells = this->schedule_
856 .changed_wells(timer.reportStepNum(), this->initialStep());
857
858 const auto changedWellLists = this->schedule_
859 .changedWellLists(timer.reportStepNum(), this->initialStep());
860
861 if (changedWells.empty() && !changedWellLists) {
862 return;
863 }
864
865 this->outputModule_->outputWellspecReport(changedWells,
866 changedWellLists,
867 timer.reportStepNum(),
868 timer.simulationTimeElapsed(),
869 timer.currentDateTime());
870 }
871
872 void writeWellflowReport(const SimulatorTimer& timer,
873 const int simStep,
874 const int wellsRequest) const
875 {
876 this->outputModule_->outputTimeStamp("WELLS",
877 timer.simulationTimeElapsed(),
878 timer.reportStepNum(),
879 timer.currentDateTime());
880
881 const auto wantConnData = wellsRequest > 1;
882
883 this->outputModule_->outputProdLog(simStep, wantConnData);
884 this->outputModule_->outputInjLog(simStep, wantConnData);
885 this->outputModule_->outputCumLog(simStep, wantConnData);
886 this->outputModule_->outputMSWLog(simStep);
887 }
888
889 int initialStep() const
890 {
891 const auto& initConfig = this->eclState().cfg().init();
892
893 return initConfig.restartRequested()
894 ? initConfig.getRestartStep()
895 : 0;
896 }
897
898 Simulator& simulator_;
899 std::unique_ptr<OutputModule> outputModule_;
900 Scalar restartTimeStepSize_;
901 int rank_ ;
902 Inplace inplace_;
903};
904
905} // namespace Opm
906
907#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: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
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 Inplace *initialInPlace, const InterRegFlowMap &interRegFlows, SummaryState &summaryState, UDQState &udqState)
Definition: EclGenericWriter_impl.hpp:1002
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