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
381 const auto start_time = boost::posix_time::
382 from_time_t(simulator_.vanguard().schedule().getStartTime());
383
384 if (this->collectOnIORank_.isIORank()) {
385 this->inplace_ = *this->outputModule_->initialInplace();
386
387 this->outputModule_->
388 outputFipAndResvLog(this->inplace_, 0, 0.0, start_time,
389 false, simulator_.gridView().comm());
390 }
391 }
392
393 outputModule_->outputFipAndResvLogToCSV(0, false, simulator_.gridView().comm());
394 }
395
396 void writeReports(const SimulatorTimer& timer)
397 {
398 if (! this->collectOnIORank_.isIORank()) {
399 return;
400 }
401
402 // SimulatorTimer::reportStepNum() is the simulator's zero-based
403 // "episode index". This is generally the index value needed to
404 // look up objects in the Schedule container. That said, function
405 // writeReports() is invoked at the *beginning* of a report
406 // step/episode which means we typically need the objects from the
407 // *previous* report step/episode. We therefore need special case
408 // handling for reportStepNum() == 0 in base runs and
409 // reportStepNum() <= restart step in restarted runs.
410 const auto firstStep = this->initialStep();
411 const auto simStep =
412 std::max(timer.reportStepNum() - 1, firstStep);
413
414 const auto& rpt = this->schedule_[simStep].rpt_config();
415
416 if (rpt.contains("WELSPECS") && (rpt.at("WELSPECS") > 0)) {
417 // Requesting a well specification report is valid at all times,
418 // including reportStepNum() == initialStep().
419 this->writeWellspecReport(timer);
420 }
421
422 if (timer.reportStepNum() == firstStep) {
423 // No dynamic flows at the beginning of the initialStep().
424 return;
425 }
426
427 if (rpt.contains("WELLS") && rpt.at("WELLS") > 0) {
428 this->writeWellflowReport(timer, simStep, rpt.at("WELLS"));
429 }
430
431 this->outputModule_->outputFipAndResvLog(this->inplace_,
432 timer.reportStepNum(),
433 timer.simulationTimeElapsed(),
434 timer.currentDateTime(),
435 /* isSubstep = */ false,
436 simulator_.gridView().comm());
437
438 OpmLog::note(""); // Blank line after all reports.
439 }
440
441 void writeOutput(data::Solution&& localCellData, bool isSubStep)
442 {
443 OPM_TIMEBLOCK(writeOutput);
444
445 const int reportStepNum = simulator_.episodeIndex() + 1;
446 this->prepareLocalCellData(isSubStep, reportStepNum);
447 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
448
449 // output using eclWriter if enabled
450 auto localWellData = simulator_.problem().wellModel().wellData();
451 auto localGroupAndNetworkData = simulator_.problem().wellModel()
452 .groupAndNetworkData(reportStepNum);
453
454 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
455 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
456
457 const bool isFlowsn = this->outputModule_->getFlows().hasFlowsn();
458 auto flowsn = this->outputModule_->getFlows().getFlowsn();
459
460 const bool isFloresn = this->outputModule_->getFlows().hasFloresn();
461 auto floresn = this->outputModule_->getFlows().getFloresn();
462
463 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
464
465 if (localCellData.empty()) {
466 this->outputModule_->assignToSolution(localCellData);
467 }
468
469 // Add cell data to perforations for RFT output
470 this->outputModule_->addRftDataToWells(localWellData,
471 reportStepNum,
472 simulator_.gridView().comm());
473 }
474
475 if (this->collectOnIORank_.isParallel() ||
476 this->collectOnIORank_.doesNeedReordering())
477 {
478 // Note: We don't need WBP (well-block averaged pressures) or
479 // inter-region flow rate values in order to create restart file
480 // output. There's consequently no need to collect those
481 // properties on the I/O rank.
482
483 this->collectOnIORank_.collect(localCellData,
484 this->outputModule_->getBlockData(),
485 this->outputModule_->getExtraBlockData(),
486 localWellData,
487 /* wbpData = */ {},
488 localGroupAndNetworkData,
489 localAquiferData,
490 localWellTestState,
491 /* interRegFlows = */ {},
492 flowsn,
493 floresn);
494 if (this->collectOnIORank_.isIORank()) {
495 this->outputModule_->assignGlobalFieldsToSolution(this->collectOnIORank_.globalCellData());
496 }
497 } else {
498 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
499 }
500
501 if (this->collectOnIORank_.isIORank()) {
502 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
503 const Scalar nextStepSize = simulator_.problem().nextTimeStepSize();
504 std::optional<int> timeStepIdx;
505 if (Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
506 timeStepIdx = simulator_.timeStepIndex();
507 }
508 this->doWriteOutput(reportStepNum, timeStepIdx, isSubStep,
509 std::move(localCellData),
510 std::move(localWellData),
511 std::move(localGroupAndNetworkData),
512 std::move(localAquiferData),
513 std::move(localWellTestState),
514 this->actionState(),
515 this->udqState(),
516 this->summaryState(),
517 this->simulator_.problem().thresholdPressure().getRestartVector(),
518 curTime, nextStepSize,
519 Parameters::Get<Parameters::EclOutputDoublePrecision>(),
520 isFlowsn, std::move(flowsn),
521 isFloresn, std::move(floresn));
522 }
523 }
524
526 {
527 const auto enablePCHysteresis = simulator_.problem().materialLawManager()->enablePCHysteresis();
528 const auto enableNonWettingHysteresis = simulator_.problem().materialLawManager()->enableNonWettingHysteresis();
529 const auto enableWettingHysteresis = simulator_.problem().materialLawManager()->enableWettingHysteresis();
530 const auto oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
531 const auto gasActive = FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
532 const auto waterActive = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
533 const auto enableSwatinit = simulator_.vanguard().eclState().fieldProps().has_double("SWATINIT");
534
535 std::vector<RestartKey> solutionKeys {
536 {"PRESSURE", UnitSystem::measure::pressure},
537 {"SWAT", UnitSystem::measure::identity, waterActive},
538 {"SGAS", UnitSystem::measure::identity, gasActive},
539 {"TEMP", UnitSystem::measure::temperature, enableEnergy},
540 {"SSOLVENT", UnitSystem::measure::identity, enableSolvent},
541
542 {"RS", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGas()},
543 {"RV", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedOil()},
544 {"RVW", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedWater()},
545 {"RSW", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGasInWater()},
546
547 {"SGMAX", UnitSystem::measure::identity, enableNonWettingHysteresis && oilActive && gasActive},
548 {"SHMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && gasActive},
549
550 {"SOMAX", UnitSystem::measure::identity,
551 (enableNonWettingHysteresis && oilActive && waterActive)
552 || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
553
554 {"SOMIN", UnitSystem::measure::identity, enablePCHysteresis && oilActive && gasActive},
555 {"SWHY1", UnitSystem::measure::identity, enablePCHysteresis && oilActive && waterActive},
556 {"SWMAX", UnitSystem::measure::identity, enableWettingHysteresis && oilActive && waterActive},
557
558 {"PPCW", UnitSystem::measure::pressure, enableSwatinit},
559 };
560
561 {
562 const auto& tracers = simulator_.vanguard().eclState().tracer();
563
564 for (const auto& tracer : tracers) {
565 const auto enableSolTracer =
566 ((tracer.phase == Phase::GAS) && FluidSystem::enableDissolvedGas()) ||
567 ((tracer.phase == Phase::OIL) && FluidSystem::enableVaporizedOil());
568
569 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity, true);
570 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
571 }
572 }
573
574 const auto& inputThpres = eclState().getSimulationConfig().getThresholdPressure();
575 const std::vector<RestartKey> extraKeys {
576 {"OPMEXTRA", UnitSystem::measure::identity, false},
577 {"THRESHPR", UnitSystem::measure::pressure, inputThpres.active()},
578 };
579
580 const auto& gridView = this->simulator_.vanguard().gridView();
581 const auto numElements = gridView.size(/*codim=*/0);
582
583 // Try to load restart step 0 to calculate initial FIP
584 {
585 this->outputModule_->allocBuffers(numElements,
586 0,
587 /*isSubStep = */false,
588 /*log = */ false,
589 /*isRestart = */true);
590
591 const auto restartSolution =
593 solutionKeys, gridView.comm(), 0);
594
595 if (!restartSolution.empty()) {
596 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
597 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
598 this->outputModule_->setRestart(restartSolution, elemIdx, globalIdx);
599 }
600
601 this->simulator_.problem().readSolutionFromOutputModule(0, true);
602 this->simulator_.problem().temperatureModel().init();
603 ElementContext elemCtx(this->simulator_);
604 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
605 elemCtx.updatePrimaryStencil(elem);
606 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
607
608 this->outputModule_->updateFluidInPlace(elemCtx);
609 }
610
611 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
612 }
613 }
614
615 {
616 // The episodeIndex is rewound one step back before calling
617 // beginRestart() and cannot be used here. We just ask the
618 // initconfig directly to be sure that we use the correct index.
619 const auto restartStepIdx = this->simulator_.vanguard()
620 .eclState().getInitConfig().getRestartStep();
621
622 this->outputModule_->allocBuffers(numElements,
623 restartStepIdx,
624 /*isSubStep = */false,
625 /*log = */ false,
626 /*isRestart = */true);
627 }
628
629 {
630 const auto restartValues =
631 loadParallelRestart(this->eclIO_.get(),
632 this->actionState(),
633 this->summaryState(),
634 solutionKeys, extraKeys, gridView.comm());
635
636 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
637 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
638 this->outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
639 }
640
641 auto& tracer_model = simulator_.problem().tracerModel();
642 for (int tracer_index = 0; tracer_index < tracer_model.numTracers(); ++tracer_index) {
643 // Free tracers
644 {
645 const auto& free_tracer_name = tracer_model.fname(tracer_index);
646 const auto& free_tracer_solution = restartValues.solution
647 .template data<double>(free_tracer_name);
648
649 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
650 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
651 tracer_model.setFreeTracerConcentration
652 (tracer_index, elemIdx, free_tracer_solution[globalIdx]);
653 }
654 }
655
656 // Solution tracer (only if DISGAS/VAPOIL are active for gas/oil tracers)
657 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
658 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil()))
659 {
660 tracer_model.setEnableSolTracers(tracer_index, true);
661
662 const auto& sol_tracer_name = tracer_model.sname(tracer_index);
663 const auto& sol_tracer_solution = restartValues.solution
664 .template data<double>(sol_tracer_name);
665
666 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
667 const auto globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
668 tracer_model.setSolTracerConcentration
669 (tracer_index, elemIdx, sol_tracer_solution[globalIdx]);
670 }
671 }
672 else {
673 tracer_model.setEnableSolTracers(tracer_index, false);
674
675 for (auto elemIdx = 0*numElements; elemIdx < numElements; ++elemIdx) {
676 tracer_model.setSolTracerConcentration(tracer_index, elemIdx, 0.0);
677 }
678 }
679 }
680
681 if (inputThpres.active()) {
682 const_cast<Simulator&>(this->simulator_)
683 .problem().thresholdPressure()
684 .setFromRestart(restartValues.getExtra("THRESHPR"));
685 }
686
687 restartTimeStepSize_ = restartValues.getExtra("OPMEXTRA")[0];
688 if (restartTimeStepSize_ <= 0) {
689 restartTimeStepSize_ = std::numeric_limits<double>::max();
690 }
691
692 // Initialize the well model from restart values
693 this->simulator_.problem().wellModel()
694 .initFromRestartFile(restartValues);
695
696 if (!restartValues.aquifer.empty()) {
697 this->simulator_.problem().mutableAquiferModel()
698 .initFromRestart(restartValues.aquifer);
699 }
700 }
701 }
702
704 {
705 // Calculate initial in-place volumes.
706 // Does nothing if they have already been calculated,
707 // e.g. from restart data at T=0.
708 this->outputModule_->calc_initial_inplace(this->simulator_.gridView().comm());
709
710 if (this->collectOnIORank_.isIORank()) {
711 if (const auto* iip = this->outputModule_->initialInplace(); iip != nullptr) {
712 this->inplace_ = *iip;
713 }
714 }
715 }
716
717 const OutputModule& outputModule() const
718 { return *outputModule_; }
719
720 OutputModule& mutableOutputModule() const
721 { return *outputModule_; }
722
723 Scalar restartTimeStepSize() const
724 { return restartTimeStepSize_; }
725
726 template <class Serializer>
727 void serializeOp(Serializer& serializer)
728 {
729 serializer(*outputModule_);
730 }
731
732private:
733 static bool enableEclOutput_()
734 {
735 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
736 return enable;
737 }
738
739 const EclipseState& eclState() const
740 { return simulator_.vanguard().eclState(); }
741
742 SummaryState& summaryState()
743 { return simulator_.vanguard().summaryState(); }
744
745 Action::State& actionState()
746 { return simulator_.vanguard().actionState(); }
747
748 UDQState& udqState()
749 { return simulator_.vanguard().udqState(); }
750
751 const Schedule& schedule() const
752 { return simulator_.vanguard().schedule(); }
753
754 void prepareLocalCellData(const bool isSubStep,
755 const int reportStepNum)
756 {
757 OPM_TIMEBLOCK(prepareLocalCellData);
758
759 if (this->outputModule_->localDataValid()) {
760 return;
761 }
762
763 const auto& gridView = simulator_.vanguard().gridView();
764 const bool log = this->collectOnIORank_.isIORank();
765
766 const int num_interior = detail::
768 this->outputModule_->
769 allocBuffers(num_interior, reportStepNum,
770 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
771 log, /*isRestart*/ false);
772
773 ElementContext elemCtx(simulator_);
774
776
777 {
778 OPM_TIMEBLOCK(prepareCellBasedData);
779
780 this->outputModule_->prepareDensityAccumulation();
781 this->outputModule_->setupExtractors(isSubStep, reportStepNum);
782 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
783 elemCtx.updatePrimaryStencil(elem);
784 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
785
786 this->outputModule_->processElement(elemCtx);
787 this->outputModule_->processElementBlockData(elemCtx);
788 }
789 this->outputModule_->clearExtractors();
790
791 this->outputModule_->accumulateDensityParallel();
792 }
793
794 {
795 OPM_TIMEBLOCK(prepareFluidInPlace);
796
797#ifdef _OPENMP
798#pragma omp parallel for
799#endif
800 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
801 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
802 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
803
804 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
805 }
806 }
807
808 this->outputModule_->validateLocalData();
809
810 OPM_END_PARALLEL_TRY_CATCH("EclWriter::prepareLocalCellData() failed: ",
811 this->simulator_.vanguard().grid().comm());
812 }
813
814 void captureLocalFluxData()
815 {
816 OPM_TIMEBLOCK(captureLocalData);
817
818 const auto& gridView = this->simulator_.vanguard().gridView();
819 const auto timeIdx = 0u;
820
821 auto elemCtx = ElementContext { this->simulator_ };
822
823 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
824 const auto activeIndex = [&elemMapper](const Element& e)
825 {
826 return elemMapper.index(e);
827 };
828
829 const auto cartesianIndex = [this](const int elemIndex)
830 {
831 return this->cartMapper_.cartesianIndex(elemIndex);
832 };
833
834 this->outputModule_->initializeFluxData();
835
837
838 for (const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
839 elemCtx.updateStencil(elem);
840 elemCtx.updateIntensiveQuantities(timeIdx);
841 elemCtx.updateExtensiveQuantities(timeIdx);
842
843 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
844 }
845
846 OPM_END_PARALLEL_TRY_CATCH("EclWriter::captureLocalFluxData() failed: ",
847 this->simulator_.vanguard().grid().comm())
848
849 this->outputModule_->finalizeFluxData();
850 }
851
852 void writeWellspecReport(const SimulatorTimer& timer) const
853 {
854 const auto changedWells = this->schedule_
855 .changed_wells(timer.reportStepNum(), this->initialStep());
856
857 const auto changedWellLists = this->schedule_
858 .changedWellLists(timer.reportStepNum(), this->initialStep());
859
860 if (changedWells.empty() && !changedWellLists) {
861 return;
862 }
863
864 this->outputModule_->outputWellspecReport(changedWells,
865 changedWellLists,
866 timer.reportStepNum(),
867 timer.simulationTimeElapsed(),
868 timer.currentDateTime());
869 }
870
871 void writeWellflowReport(const SimulatorTimer& timer,
872 const int simStep,
873 const int wellsRequest) const
874 {
875 this->outputModule_->outputTimeStamp("WELLS",
876 timer.simulationTimeElapsed(),
877 timer.reportStepNum(),
878 timer.currentDateTime());
879
880 const auto wantConnData = wellsRequest > 1;
881
882 this->outputModule_->outputProdLog(simStep, wantConnData);
883 this->outputModule_->outputInjLog(simStep, wantConnData);
884 this->outputModule_->outputCumLog(simStep, wantConnData);
885 this->outputModule_->outputMSWLog(simStep);
886 }
887
888 int initialStep() const
889 {
890 const auto& initConfig = this->eclState().cfg().init();
891
892 return initConfig.restartRequested()
893 ? initConfig.getRestartStep()
894 : 0;
895 }
896
897 Simulator& simulator_;
898 std::unique_ptr<OutputModule> outputModule_;
899 Scalar restartTimeStepSize_;
900 int rank_ ;
901 Inplace inplace_;
902};
903
904} // namespace Opm
905
906#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:66
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:553
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:661
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:115
OutputModule & mutableOutputModule() const
Definition: EclWriter.hpp:720
const OutputModule & outputModule() const
Definition: EclWriter.hpp:717
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:727
void writeInitialFIPReport()
Writes the initial FIP report as configured in RPTSOL.
Definition: EclWriter.hpp:351
void beginRestart()
Definition: EclWriter.hpp:525
EclWriter(Simulator &simulator)
Definition: EclWriter.hpp:152
void writeOutput(data::Solution &&localCellData, bool isSubStep)
Definition: EclWriter.hpp:441
void writeReports(const SimulatorTimer &timer)
Definition: EclWriter.hpp:396
void endRestart()
Definition: EclWriter.hpp:703
Scalar restartTimeStepSize() const
Definition: EclWriter.hpp:723
~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