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
57#ifdef RESERVOIR_COUPLING_ENABLED
59#endif
60
61#include <boost/date_time/posix_time/posix_time.hpp>
62
63#include <limits>
64#include <map>
65#include <memory>
66#include <optional>
67#include <stdexcept>
68#include <string>
69#include <utility>
70#include <vector>
71
72namespace Opm::Parameters {
73
74// If available, write the ECL output in a non-blocking manner
75struct EnableAsyncEclOutput { static constexpr bool value = true; };
76
77// By default, use single precision for the ECL formated results
78struct EclOutputDoublePrecision { static constexpr bool value = false; };
79
80// Write all solutions for visualization, not just the ones for the
81// report steps...
82struct EnableWriteAllSolutions { static constexpr bool value = false; };
83
84// Write ESMRY file for fast loading of summary data
85struct EnableEsmry { static constexpr bool value = true; };
86
87} // namespace Opm::Parameters
88
89namespace Opm::Action {
90 class State;
91} // namespace Opm::Action
92
93namespace Opm {
94 class EclipseIO;
95 class UDQState;
96} // namespace Opm
97
98namespace Opm {
114template <class TypeTag, class OutputModule>
115class EclWriter : public EclGenericWriter<GetPropType<TypeTag, Properties::Grid>,
116 GetPropType<TypeTag, Properties::EquilGrid>,
117 GetPropType<TypeTag, Properties::GridView>,
118 GetPropType<TypeTag, Properties::ElementMapper>,
119 GetPropType<TypeTag, Properties::Scalar>>
120{
129 using Element = typename GridView::template Codim<0>::Entity;
131 using ElementIterator = typename GridView::template Codim<0>::Iterator;
133
134 typedef Dune::MultipleCodimMultipleGeomTypeMapper< GridView > VertexMapper;
135
136 enum { enableEnergy = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal ||
137 getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::SequentialImplicitThermal };
138 enum { enableMech = getPropValue<TypeTag, Properties::EnableMech>() };
139 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
140 enum { enableGeochemistry = getPropValue<TypeTag, Properties::EnableGeochemistry>() };
141
142public:
143
144 static void registerParameters()
145 {
146 OutputModule::registerParameters();
147
148 Parameters::Register<Parameters::EnableAsyncEclOutput>
149 ("Write the ECL-formated results in a non-blocking way "
150 "(i.e., using a separate thread).");
151 Parameters::Register<Parameters::EnableEsmry>
152 ("Write ESMRY file for fast loading of summary data.");
153 }
154
155 // The Simulator object should preferably have been const - the
156 // only reason that is not the case is due to the SummaryState
157 // object owned deep down by the vanguard.
158 explicit EclWriter(Simulator& simulator)
159 : BaseType(simulator.vanguard().schedule(),
160 simulator.vanguard().eclState(),
161 simulator.vanguard().summaryConfig(),
162 simulator.vanguard().grid(),
163 ((simulator.vanguard().grid().comm().rank() == 0)
164 ? &simulator.vanguard().equilGrid()
165 : nullptr),
166 simulator.vanguard().gridView(),
167 simulator.vanguard().cartesianIndexMapper(),
168 ((simulator.vanguard().grid().comm().rank() == 0)
169 ? &simulator.vanguard().equilCartesianIndexMapper()
170 : nullptr),
171 Parameters::Get<Parameters::EnableAsyncEclOutput>(),
172 Parameters::Get<Parameters::EnableEsmry>())
173 , simulator_(simulator)
174 {
175#if HAVE_MPI
176 if (this->simulator_.vanguard().grid().comm().size() > 1) {
177 auto smryCfg = (this->simulator_.vanguard().grid().comm().rank() == 0)
178 ? this->eclIO_->finalSummaryConfig()
179 : SummaryConfig{};
180
181 eclBroadcast(this->simulator_.vanguard().grid().comm(), smryCfg);
182
183 this->outputModule_ = std::make_unique<OutputModule>
184 (simulator, smryCfg, this->collectOnIORank_);
185 }
186 else
187#endif
188 {
189 this->outputModule_ = std::make_unique<OutputModule>
190 (simulator, this->eclIO_->finalSummaryConfig(), this->collectOnIORank_);
191 }
192
193 this->rank_ = this->simulator_.vanguard().grid().comm().rank();
194
195 this->simulator_.vanguard().eclState().computeFipRegionStatistics();
196 }
197
199 {}
200
201 const EquilGrid& globalGrid() const
202 {
203 return simulator_.vanguard().equilGrid();
204 }
205
209 void evalSummaryState(bool isSubStep)
210 {
211 OPM_TIMEBLOCK(evalSummaryState);
212 const int reportStepNum = simulator_.episodeIndex() + 1;
213
214 /*
215 The summary data is not evaluated for timestep 0, that is
216 implemented with a:
217
218 if (time_step == 0)
219 return;
220
221 check somewhere in the summary code. When the summary code was
222 split in separate methods Summary::eval() and
223 Summary::add_timestep() it was necessary to pull this test out
224 here to ensure that the well and group related keywords in the
225 restart file, like XWEL and XGRP were "correct" also in the
226 initial report step.
227
228 "Correct" in this context means unchanged behavior, might very
229 well be more correct to actually remove this if test.
230 */
231
232 if (reportStepNum == 0)
233 return;
234
235 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
236 const Scalar totalCpuTime =
237 simulator_.executionTimer().realTimeElapsed() +
238 simulator_.setupTimer().realTimeElapsed() +
239 simulator_.vanguard().setupTime();
240
241 const auto localWellData = simulator_.problem().wellModel().wellData();
242 const auto localWBP = simulator_.problem().wellModel().wellBlockAveragePressures();
243 const auto localGroupAndNetworkData = simulator_.problem().wellModel()
244 .groupAndNetworkData(reportStepNum);
245
246 const auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
247 const auto localWellTestState = simulator_.problem().wellModel().wellTestState();
248 this->prepareLocalCellData(isSubStep, reportStepNum);
249
250 if (this->outputModule_->needInterfaceFluxes(isSubStep)) {
251 this->captureLocalFluxData();
252 }
253
254 if (this->collectOnIORank_.isParallel()) {
256
257 std::map<std::pair<std::string,int>,double> dummy;
258 this->collectOnIORank_.collect({},
259 outputModule_->getBlockData(),
260 dummy,
261 localWellData,
262 localWBP,
263 localGroupAndNetworkData,
264 localAquiferData,
265 localWellTestState,
266 this->outputModule_->getInterRegFlows(),
267 {},
268 {});
269
270 if (this->collectOnIORank_.isIORank()) {
271 auto& iregFlows = this->collectOnIORank_.globalInterRegFlows();
272
273 if (! iregFlows.readIsConsistent()) {
274 throw std::runtime_error {
275 "Inconsistent inter-region flow "
276 "region set names in parallel"
277 };
278 }
279
280 iregFlows.compress();
281 }
282
283 OPM_END_PARALLEL_TRY_CATCH("Collect to I/O rank: ",
284 this->simulator_.vanguard().grid().comm());
285 }
286
287
288 std::map<std::string, double> miscSummaryData;
289 std::map<std::string, std::vector<double>> regionData;
290 Inplace inplace;
291
292 {
293 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
294
295 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
296
297 if (this->collectOnIORank_.isIORank()){
298 inplace_ = inplace;
299 }
300 }
301
302 // Add TCPU
303 if (totalCpuTime != 0.0) {
304 miscSummaryData["TCPU"] = totalCpuTime;
305 }
307 miscSummaryData["NEWTON"] = this->sub_step_report_.total_newton_iterations;
308 }
310 miscSummaryData["MLINEARS"] = this->sub_step_report_.total_linear_iterations;
311 }
313 miscSummaryData["NLINEARS"] = static_cast<float>(this->sub_step_report_.total_linear_iterations) / this->sub_step_report_.total_newton_iterations;
314 }
315 if (this->sub_step_report_.min_linear_iterations != std::numeric_limits<unsigned int>::max()) {
316 miscSummaryData["NLINSMIN"] = this->sub_step_report_.min_linear_iterations;
317 }
319 miscSummaryData["NLINSMAX"] = this->sub_step_report_.max_linear_iterations;
320 }
322 miscSummaryData["MSUMLINS"] = this->simulation_report_.success.total_linear_iterations;
323 }
325 miscSummaryData["MSUMNEWT"] = this->simulation_report_.success.total_newton_iterations;
326 }
327
328 // For reservoir coupling master: collect slave production/injection
329 // rates to pass through to Summary::eval() via DynamicSimulatorState.
330 const auto rcGroupRates = this->collectReservoirCouplingGroupRates_();
331
332 {
333 OPM_TIMEBLOCK(evalSummary);
334
335 const auto& blockData = this->collectOnIORank_.isParallel()
337 : this->outputModule_->getBlockData();
338
339 const auto& interRegFlows = this->collectOnIORank_.isParallel()
341 : this->outputModule_->getInterRegFlows();
342
343 this->evalSummary(reportStepNum,
344 curTime,
345 localWellData,
346 localWBP,
347 localGroupAndNetworkData,
348 localAquiferData,
349 blockData,
350 miscSummaryData,
351 regionData,
352 inplace,
353 this->outputModule_->initialInplace(),
354 interRegFlows,
355 this->summaryState(),
356 this->udqState(),
357 rcGroupRates ? &(*rcGroupRates) : nullptr);
358 }
359 }
360
363 {
364 const auto& gridView = simulator_.vanguard().gridView();
365 const int num_interior = detail::
367
368 this->outputModule_->
369 allocBuffers(num_interior, 0, false, false, /*isRestart*/ false);
370
371#ifdef _OPENMP
372#pragma omp parallel for
373#endif
374 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
375 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
376 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
377
378 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
379 }
380
381 // We always calculate the initial fip values as it may be used by various
382 // keywords in the Schedule, e.g. FIP=2 in RPTSCHED but no FIP in RPTSOL
383 outputModule_->calc_initial_inplace(simulator_.gridView().comm());
384
385 // check if RPTSOL entry has FIP output
386 const auto& fip = simulator_.vanguard().eclState().getEclipseConfig().fip();
387 if (fip.output(FIPConfig::OutputField::FIELD) ||
388 fip.output(FIPConfig::OutputField::RESV))
389 {
390 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
391
392 const auto start_time = boost::posix_time::
393 from_time_t(simulator_.vanguard().schedule().getStartTime());
394
395 if (this->collectOnIORank_.isIORank()) {
396 this->inplace_ = *this->outputModule_->initialInplace();
397
398 this->outputModule_->
399 outputFipAndResvLog(this->inplace_, 0, 0.0, start_time,
400 false, simulator_.gridView().comm());
401 }
402 }
403
404 outputModule_->outputFipAndResvLogToCSV(0, false, simulator_.gridView().comm());
405 }
406
407 void writeReports(const SimulatorTimer& timer)
408 {
409 if (! this->collectOnIORank_.isIORank()) {
410 return;
411 }
412
413 // SimulatorTimer::reportStepNum() is the simulator's zero-based
414 // "episode index". This is generally the index value needed to
415 // look up objects in the Schedule container. That said, function
416 // writeReports() is invoked at the *beginning* of a report
417 // step/episode which means we typically need the objects from the
418 // *previous* report step/episode. We therefore need special case
419 // handling for reportStepNum() == 0 in base runs and
420 // reportStepNum() <= restart step in restarted runs.
421 const auto firstStep = this->initialStep();
422 const auto simStep =
423 std::max(timer.reportStepNum() - 1, firstStep);
424
425 const auto& rpt = this->schedule_[simStep].rpt_config();
426
427 if (rpt.contains("WELSPECS") && (rpt.at("WELSPECS") > 0)) {
428 // Requesting a well specification report is valid at all times,
429 // including reportStepNum() == initialStep().
430 this->writeWellspecReport(timer);
431 }
432
433 if (timer.reportStepNum() == firstStep) {
434 // No dynamic flows at the beginning of the initialStep().
435 return;
436 }
437
438 if (rpt.contains("WELLS") && rpt.at("WELLS") > 0) {
439 this->writeWellflowReport(timer, simStep, rpt.at("WELLS"));
440 }
441
442 this->outputModule_->outputFipAndResvLog(this->inplace_,
443 timer.reportStepNum(),
444 timer.simulationTimeElapsed(),
445 timer.currentDateTime(),
446 /* isSubstep = */ false,
447 simulator_.gridView().comm());
448
449 OpmLog::note(""); // Blank line after all reports.
450 }
451
452 void writeOutput(data::Solution&& localCellData, const bool isSubStep, const bool isForcedFinalOutput)
453 {
454 OPM_TIMEBLOCK(writeOutput);
455
456 const int reportStepNum = simulator_.episodeIndex() + 1;
457 this->prepareLocalCellData(isSubStep, reportStepNum);
458 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
459
460 // output using eclWriter if enabled
461 auto localWellData = simulator_.problem().wellModel().wellData();
462 auto localGroupAndNetworkData = simulator_.problem().wellModel()
463 .groupAndNetworkData(reportStepNum);
464
465 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
466 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
467
468 const bool isFlowsn = this->outputModule_->getFlows().hasFlowsn();
469 auto flowsn = this->outputModule_->getFlows().getFlowsn();
470
471 const bool isFloresn = this->outputModule_->getFlows().hasFloresn();
472 auto floresn = this->outputModule_->getFlows().getFloresn();
473
474 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
475
476 if (localCellData.empty()) {
477 this->outputModule_->assignToSolution(localCellData);
478 }
479
480 // Add cell data to perforations for RFT output
481 this->outputModule_->addRftDataToWells(localWellData,
482 reportStepNum,
483 simulator_.gridView().comm());
484 }
485
486 if (this->collectOnIORank_.isParallel() ||
487 this->collectOnIORank_.doesNeedReordering())
488 {
489 // Note: We don't need WBP (well-block averaged pressures) or
490 // inter-region flow rate values in order to create restart file
491 // output. There's consequently no need to collect those
492 // properties on the I/O rank.
493
494 this->collectOnIORank_.collect(localCellData,
495 this->outputModule_->getBlockData(),
496 this->outputModule_->getExtraBlockData(),
497 localWellData,
498 /* wbpData = */ {},
499 localGroupAndNetworkData,
500 localAquiferData,
501 localWellTestState,
502 /* interRegFlows = */ {},
503 flowsn,
504 floresn);
505 if (this->collectOnIORank_.isIORank()) {
506 this->outputModule_->assignGlobalFieldsToSolution(this->collectOnIORank_.globalCellData());
507 }
508 } else {
509 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
510 }
511
512 if (this->collectOnIORank_.isIORank()) {
513 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
514 const Scalar nextStepSize = simulator_.problem().nextTimeStepSize();
515 std::optional<int> timeStepIdx;
516 if (Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
517 timeStepIdx = simulator_.timeStepIndex();
518 }
519 this->doWriteOutput(reportStepNum, timeStepIdx, isSubStep,
520 isForcedFinalOutput,
521 std::move(localCellData),
522 std::move(localWellData),
523 std::move(localGroupAndNetworkData),
524 std::move(localAquiferData),
525 std::move(localWellTestState),
526 this->actionState(),
527 this->udqState(),
528 this->summaryState(),
529 this->simulator_.problem().thresholdPressure().getRestartVector(),
530 curTime, nextStepSize,
531 Parameters::Get<Parameters::EclOutputDoublePrecision>(),
532 isFlowsn, std::move(flowsn),
533 isFloresn, std::move(floresn));
534 }
535 }
536
538 {
539 const auto enablePCHysteresis = simulator_.problem().materialLawManager()->enablePCHysteresis();
540 const auto enableNonWettingHysteresis = simulator_.problem().materialLawManager()->enableNonWettingHysteresis();
541 const auto enableWettingHysteresis = simulator_.problem().materialLawManager()->enableWettingHysteresis();
542 const auto oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
543 const auto gasActive = FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
544 const auto waterActive = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
545 const auto enableSwatinit = simulator_.vanguard().eclState().fieldProps().has_double("SWATINIT");
546
547 std::vector<RestartKey> solutionKeys {
548 {"PRESSURE", UnitSystem::measure::pressure},
549 {"SWAT", UnitSystem::measure::identity, waterActive},
550 {"SGAS", UnitSystem::measure::identity, gasActive},
551 {"TEMP", UnitSystem::measure::temperature, enableEnergy},
552 {"SSOLVENT", UnitSystem::measure::identity, enableSolvent},
553
554 {"RS", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGas()},
555 {"RV", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedOil()},
556 {"RVW", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedWater()},
557 {"RSW", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGasInWater()},
558
559 {"SGMAX", UnitSystem::measure::identity, enableNonWettingHysteresis && oilActive && gasActive},
560 {"SHMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && gasActive},
561
562 {"SOMAX", UnitSystem::measure::identity,
563 (enableNonWettingHysteresis && oilActive && waterActive)
564 || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
565
566 {"SOMIN", UnitSystem::measure::identity, enablePCHysteresis && oilActive && gasActive},
567 {"SWHY1", UnitSystem::measure::identity, enablePCHysteresis && oilActive && waterActive},
568 {"SWMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && waterActive},
569
570 {"PPCW", UnitSystem::measure::pressure, enableSwatinit},
571 };
572
573 {
574 const auto& tracers = simulator_.vanguard().eclState().tracer();
575
576 for (const auto& tracer : tracers) {
577 const auto enableSolTracer =
578 ((tracer.phase == Phase::GAS) && FluidSystem::enableDissolvedGas()) ||
579 ((tracer.phase == Phase::OIL) && FluidSystem::enableVaporizedOil());
580
581 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity, true);
582 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
583 }
584 }
585
586 const auto& inputThpres = eclState().getSimulationConfig().getThresholdPressure();
587 const std::vector<RestartKey> extraKeys {
588 {"OPMEXTRA", UnitSystem::measure::identity, false},
589 {"THRESHPR", UnitSystem::measure::pressure, inputThpres.active()},
590 };
591
592 const auto& gridView = this->simulator_.vanguard().gridView();
593 const auto numElements = gridView.size(/*codim=*/0);
594
595 // Try to load restart step 0 to calculate initial FIP
596 {
597 this->outputModule_->allocBuffers(numElements,
598 0,
599 /*isSubStep = */false,
600 /*log = */ false,
601 /*isRestart = */true);
602
603 const auto restartSolution =
605 solutionKeys, gridView.comm(), 0);
606
607 if (!restartSolution.empty()) {
608 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
609 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
610 this->outputModule_->setRestart(restartSolution, elemIdx, globalIdx);
611 }
612
613 this->simulator_.problem().readSolutionFromOutputModule(0, true);
614 this->simulator_.problem().temperatureModel().init();
615 ElementContext elemCtx(this->simulator_);
616 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
617 elemCtx.updatePrimaryStencil(elem);
618 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
619
620 this->outputModule_->updateFluidInPlace(elemCtx);
621 }
622
623 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
624 }
625 }
626
627 {
628 // The episodeIndex is rewound one step back before calling
629 // beginRestart() and cannot be used here. We just ask the
630 // initconfig directly to be sure that we use the correct index.
631 const auto restartStepIdx = this->simulator_.vanguard()
632 .eclState().getInitConfig().getRestartStep();
633
634 this->outputModule_->allocBuffers(numElements,
635 restartStepIdx,
636 /*isSubStep = */false,
637 /*log = */ false,
638 /*isRestart = */true);
639 }
640
641 {
642 const auto restartValues =
643 loadParallelRestart(this->eclIO_.get(),
644 this->actionState(),
645 this->summaryState(),
646 solutionKeys, extraKeys, gridView.comm());
647
648 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
649 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
650 this->outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
651 }
652
653 auto& tracer_model = simulator_.problem().tracerModel();
654 for (int tracer_index = 0; tracer_index < tracer_model.numTracers(); ++tracer_index) {
655 // Free tracers
656 {
657 const auto& free_tracer_name = tracer_model.fname(tracer_index);
658 const auto& free_tracer_solution = restartValues.solution
659 .template data<double>(free_tracer_name);
660
661 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
662 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
663 tracer_model.setFreeTracerConcentration
664 (tracer_index, elemIdx, free_tracer_solution[globalIdx]);
665 }
666 }
667
668 // Solution tracer (only if DISGAS/VAPOIL are active for gas/oil tracers)
669 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
670 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil()))
671 {
672 tracer_model.setEnableSolTracers(tracer_index, true);
673
674 const auto& sol_tracer_name = tracer_model.sname(tracer_index);
675 const auto& sol_tracer_solution = restartValues.solution
676 .template data<double>(sol_tracer_name);
677
678 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
679 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
680 tracer_model.setSolTracerConcentration
681 (tracer_index, elemIdx, sol_tracer_solution[globalIdx]);
682 }
683 }
684 else {
685 tracer_model.setEnableSolTracers(tracer_index, false);
686
687 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
688 tracer_model.setSolTracerConcentration(tracer_index, elemIdx, 0.0);
689 }
690 }
691 }
692
693 if (inputThpres.active()) {
694 const_cast<Simulator&>(this->simulator_)
695 .problem().thresholdPressure()
696 .setFromRestart(restartValues.getExtra("THRESHPR"));
697 }
698
699 restartTimeStepSize_ = restartValues.getExtra("OPMEXTRA")[0];
700 if (restartTimeStepSize_ <= 0) {
701 restartTimeStepSize_ = std::numeric_limits<double>::max();
702 }
703
704 // Initialize the well model from restart values
705 this->simulator_.problem().wellModel()
706 .initFromRestartFile(restartValues);
707
708 if (!restartValues.aquifer.empty()) {
709 this->simulator_.problem().mutableAquiferModel()
710 .initFromRestart(restartValues.aquifer);
711 }
712 }
713 }
714
716 {
717 // Calculate initial in-place volumes.
718 // Does nothing if they have already been calculated,
719 // e.g. from restart data at T=0.
720 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
721
722 if (this->collectOnIORank_.isIORank()) {
723 if (const auto* iip = this->outputModule_->initialInplace(); iip != nullptr) {
724 this->inplace_ = *iip;
725 }
726 }
727 }
728
729 const OutputModule& outputModule() const
730 { return *outputModule_; }
731
732 OutputModule& mutableOutputModule() const
733 { return *outputModule_; }
734
735 Scalar restartTimeStepSize() const
736 { return restartTimeStepSize_; }
737
738 template <class Serializer>
739 void serializeOp(Serializer& serializer)
740 {
741 serializer(*outputModule_);
742 }
743
744private:
745 static bool enableEclOutput_()
746 {
747 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
748 return enable;
749 }
750
751 const EclipseState& eclState() const
752 { return simulator_.vanguard().eclState(); }
753
754 SummaryState& summaryState()
755 { return simulator_.vanguard().summaryState(); }
756
757 Action::State& actionState()
758 { return simulator_.vanguard().actionState(); }
759
760 UDQState& udqState()
761 { return simulator_.vanguard().udqState(); }
762
763 const Schedule& schedule() const
764 { return simulator_.vanguard().schedule(); }
765
768 std::optional<data::ReservoirCouplingGroupRates> collectReservoirCouplingGroupRates_()
769 {
770#ifdef RESERVOIR_COUPLING_ENABLED
771 // Guard: only BlackoilWellModel has reservoir coupling support.
772 // CompWellModel (compositional) does not, so we use if constexpr
773 // to avoid compilation errors when EclWriter is instantiated with
774 // a compositional TypeTag.
775 using WellModelType = std::remove_cvref_t<
776 decltype(simulator_.problem().wellModel())>;
777 if constexpr (requires(WellModelType& wm) { wm.isReservoirCouplingMaster(); }) {
778 auto& wellModel = simulator_.problem().wellModel();
779 if (!wellModel.isReservoirCouplingMaster()) {
780 return std::nullopt;
781 }
782 return wellModel.reservoirCouplingMaster()
783 .collectGroupRatesForSummary();
784 }
785#endif
786 return std::nullopt;
787 }
788
789 void prepareLocalCellData(const bool isSubStep,
790 const int reportStepNum)
791 {
792 OPM_TIMEBLOCK(prepareLocalCellData);
793
794 if (this->outputModule_->localDataValid()) {
795 return;
796 }
797
798 const auto& gridView = simulator_.vanguard().gridView();
799 const bool log = this->collectOnIORank_.isIORank();
800
801 const int num_interior = detail::
803 this->outputModule_->
804 allocBuffers(num_interior, reportStepNum,
805 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
806 log, /*isRestart*/ false);
807
808 ElementContext elemCtx(simulator_);
809
811
812 {
813 OPM_TIMEBLOCK(prepareCellBasedData);
814
815 this->outputModule_->prepareDensityAccumulation();
816 this->outputModule_->setupExtractors(isSubStep, reportStepNum);
817 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
818 elemCtx.updatePrimaryStencil(elem);
819 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
820
821 this->outputModule_->processElement(elemCtx);
822 this->outputModule_->processElementBlockData(elemCtx);
823 }
824 this->outputModule_->clearExtractors();
825
826 this->outputModule_->accumulateDensityParallel();
827 }
828
829 {
830 OPM_TIMEBLOCK(prepareFluidInPlace);
831
832#ifdef _OPENMP
833#pragma omp parallel for
834#endif
835 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
836 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
837 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
838
839 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
840 }
841 }
842
843 this->outputModule_->validateLocalData();
844
845 OPM_END_PARALLEL_TRY_CATCH("EclWriter::prepareLocalCellData() failed: ",
846 this->simulator_.vanguard().grid().comm());
847 }
848
849 void captureLocalFluxData()
850 {
851 OPM_TIMEBLOCK(captureLocalData);
852
853 const auto& gridView = this->simulator_.vanguard().gridView();
854 const auto timeIdx = 0u;
855
856 auto elemCtx = ElementContext { this->simulator_ };
857
858 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
859 const auto activeIndex = [&elemMapper](const Element& e)
860 {
861 return elemMapper.index(e);
862 };
863
864 const auto cartesianIndex = [this](const int elemIndex)
865 {
866 return this->cartMapper_.cartesianIndex(elemIndex);
867 };
868
869 this->outputModule_->initializeFluxData();
870
872
873 for (const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
874 elemCtx.updateStencil(elem);
875 elemCtx.updateIntensiveQuantities(timeIdx);
876 elemCtx.updateExtensiveQuantities(timeIdx);
877
878 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
879 }
880
881 OPM_END_PARALLEL_TRY_CATCH("EclWriter::captureLocalFluxData() failed: ",
882 this->simulator_.vanguard().grid().comm())
883
884 this->outputModule_->finalizeFluxData();
885 }
886
887 void writeWellspecReport(const SimulatorTimer& timer) const
888 {
889 const auto changedWells = this->schedule_
890 .changed_wells(timer.reportStepNum(), this->initialStep());
891
892 const auto changedWellLists = this->schedule_
893 .changedWellLists(timer.reportStepNum(), this->initialStep());
894
895 if (changedWells.empty() && !changedWellLists) {
896 return;
897 }
898
899 this->outputModule_->outputWellspecReport(changedWells,
900 changedWellLists,
901 timer.reportStepNum(),
902 timer.simulationTimeElapsed(),
903 timer.currentDateTime());
904 }
905
906 void writeWellflowReport(const SimulatorTimer& timer,
907 const int simStep,
908 const int wellsRequest) const
909 {
910 this->outputModule_->outputTimeStamp("WELLS",
911 timer.simulationTimeElapsed(),
912 timer.reportStepNum(),
913 timer.currentDateTime());
914
915 const auto wantConnData = wellsRequest > 1;
916
917 this->outputModule_->outputProdLog(simStep, wantConnData);
918 this->outputModule_->outputInjLog(simStep, wantConnData);
919 this->outputModule_->outputCumLog(simStep, wantConnData);
920 this->outputModule_->outputMSWLog(simStep);
921 }
922
923 int initialStep() const
924 {
925 const auto& initConfig = this->eclState().cfg().init();
926
927 return initConfig.restartRequested()
928 ? initConfig.getRestartStep()
929 : 0;
930 }
931
932 Simulator& simulator_;
933 std::unique_ptr<OutputModule> outputModule_;
934 Scalar restartTimeStepSize_;
935 int rank_ ;
936 Inplace inplace_;
937};
938
939} // namespace Opm
940
941#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:69
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, const data::ReservoirCouplingGroupRates *rcGroupRates=nullptr)
Definition: EclGenericWriter_impl.hpp:1021
void doWriteOutput(const int reportStepNum, const std::optional< int > timeStepNum, const bool isSubStep, const bool forcedSimulationFinished, 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:911
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:120
OutputModule & mutableOutputModule() const
Definition: EclWriter.hpp:732
const OutputModule & outputModule() const
Definition: EclWriter.hpp:729
void writeOutput(data::Solution &&localCellData, const bool isSubStep, const bool isForcedFinalOutput)
Definition: EclWriter.hpp:452
void evalSummaryState(bool isSubStep)
collect and pass data and pass it to eclIO writer
Definition: EclWriter.hpp:209
static void registerParameters()
Definition: EclWriter.hpp:144
void serializeOp(Serializer &serializer)
Definition: EclWriter.hpp:739
void writeInitialFIPReport()
Writes the initial FIP report as configured in RPTSOL.
Definition: EclWriter.hpp:362
void beginRestart()
Definition: EclWriter.hpp:537
EclWriter(Simulator &simulator)
Definition: EclWriter.hpp:158
void writeReports(const SimulatorTimer &timer)
Definition: EclWriter.hpp:407
void endRestart()
Definition: EclWriter.hpp:715
Scalar restartTimeStepSize() const
Definition: EclWriter.hpp:735
~EclWriter()
Definition: EclWriter.hpp:198
const EquilGrid & globalGrid() const
Definition: EclWriter.hpp:201
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:78
static constexpr bool value
Definition: EclWriter.hpp:78
Definition: EclWriter.hpp:75
static constexpr bool value
Definition: EclWriter.hpp:75
Definition: EclWriter.hpp:85
static constexpr bool value
Definition: EclWriter.hpp:85
Definition: EclWriter.hpp:82
static constexpr bool value
Definition: EclWriter.hpp:82
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