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/OpmLog/OpmLog.hpp>
34#include <opm/input/eclipse/Schedule/RPTConfig.hpp>
35
36#include <opm/input/eclipse/Units/UnitSystem.hpp>
37
38#include <opm/output/eclipse/RestartValue.hpp>
39
49
50#include <limits>
51#include <stdexcept>
52#include <string>
53
54namespace Opm::Parameters {
55
56// enable the ECL output by default
57struct EnableEclOutput { static constexpr bool value = true; };
58
59// If available, write the ECL output in a non-blocking manner
60struct EnableAsyncEclOutput { static constexpr bool value = true; };
61
62// By default, use single precision for the ECL formated results
63struct EclOutputDoublePrecision { static constexpr bool value = false; };
64
65// Write all solutions for visualization, not just the ones for the
66// report steps...
67struct EnableWriteAllSolutions { static constexpr bool value = false; };
68
69// Write ESMRY file for fast loading of summary data
70struct EnableEsmry { static constexpr bool value = false; };
71
72} // namespace Opm::Parameters
73
74namespace Opm {
75
76namespace Action { class State; }
77class EclipseIO;
78class UDQState;
79
95template <class TypeTag>
96class EclWriter : public EclGenericWriter<GetPropType<TypeTag, Properties::Grid>,
97 GetPropType<TypeTag, Properties::EquilGrid>,
98 GetPropType<TypeTag, Properties::GridView>,
99 GetPropType<TypeTag, Properties::ElementMapper>,
100 GetPropType<TypeTag, Properties::Scalar>>
101{
110 using Element = typename GridView::template Codim<0>::Entity;
112 using ElementIterator = typename GridView::template Codim<0>::Iterator;
114
115 typedef Dune::MultipleCodimMultipleGeomTypeMapper< GridView > VertexMapper;
116
117 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
118 enum { enableMech = getPropValue<TypeTag, Properties::EnableMech>() };
119 enum { enableTemperature = getPropValue<TypeTag, Properties::EnableTemperature>() };
120 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
121
122public:
123
124 static void registerParameters()
125 {
127
128 Parameters::Register<Parameters::EnableAsyncEclOutput>
129 ("Write the ECL-formated results in a non-blocking way "
130 "(i.e., using a separate thread).");
131 Parameters::Register<Parameters::EnableEsmry>
132 ("Write ESMRY file for fast loading of summary data.");
133 }
134
135 // The Simulator object should preferably have been const - the
136 // only reason that is not the case is due to the SummaryState
137 // object owned deep down by the vanguard.
138 EclWriter(Simulator& simulator)
139 : BaseType(simulator.vanguard().schedule(),
140 simulator.vanguard().eclState(),
141 simulator.vanguard().summaryConfig(),
142 simulator.vanguard().grid(),
143 ((simulator.vanguard().grid().comm().rank() == 0)
144 ? &simulator.vanguard().equilGrid()
145 : nullptr),
146 simulator.vanguard().gridView(),
147 simulator.vanguard().cartesianIndexMapper(),
148 ((simulator.vanguard().grid().comm().rank() == 0)
149 ? &simulator.vanguard().equilCartesianIndexMapper()
150 : nullptr),
151 Parameters::Get<Parameters::EnableAsyncEclOutput>(),
152 Parameters::Get<Parameters::EnableEsmry>())
153 , simulator_(simulator)
154 {
155#if HAVE_MPI
156 if (this->simulator_.vanguard().grid().comm().size() > 1) {
157 auto smryCfg = (this->simulator_.vanguard().grid().comm().rank() == 0)
158 ? this->eclIO_->finalSummaryConfig()
159 : SummaryConfig{};
160
161 eclBroadcast(this->simulator_.vanguard().grid().comm(), smryCfg);
162
163 this->outputModule_ = std::make_unique<OutputBlackOilModule<TypeTag>>
164 (simulator, smryCfg, this->collectOnIORank_);
165 }
166 else
167#endif
168 {
169 this->outputModule_ = std::make_unique<OutputBlackOilModule<TypeTag>>
170 (simulator, this->eclIO_->finalSummaryConfig(), this->collectOnIORank_);
171 }
172
173 this->rank_ = this->simulator_.vanguard().grid().comm().rank();
174
175 this->simulator_.vanguard().eclState().computeFipRegionStatistics();
176 }
177
179 {}
180
181 const EquilGrid& globalGrid() const
182 {
183 return simulator_.vanguard().equilGrid();
184 }
185
189 void evalSummaryState(bool isSubStep)
190 {
191 OPM_TIMEBLOCK(evalSummaryState);
192 const int reportStepNum = simulator_.episodeIndex() + 1;
193
194 /*
195 The summary data is not evaluated for timestep 0, that is
196 implemented with a:
197
198 if (time_step == 0)
199 return;
200
201 check somewhere in the summary code. When the summary code was
202 split in separate methods Summary::eval() and
203 Summary::add_timestep() it was necessary to pull this test out
204 here to ensure that the well and group related keywords in the
205 restart file, like XWEL and XGRP were "correct" also in the
206 initial report step.
207
208 "Correct" in this context means unchanged behavior, might very
209 well be more correct to actually remove this if test.
210 */
211
212 if (reportStepNum == 0)
213 return;
214
215 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
216 const Scalar totalCpuTime =
217 simulator_.executionTimer().realTimeElapsed() +
218 simulator_.setupTimer().realTimeElapsed() +
219 simulator_.vanguard().setupTime();
220
221 const auto localWellData = simulator_.problem().wellModel().wellData();
222 const auto localWBP = simulator_.problem().wellModel().wellBlockAveragePressures();
223 const auto localGroupAndNetworkData = simulator_.problem().wellModel()
224 .groupAndNetworkData(reportStepNum);
225
226 const auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
227 const auto localWellTestState = simulator_.problem().wellModel().wellTestState();
228 this->prepareLocalCellData(isSubStep, reportStepNum);
229
230 if (this->outputModule_->needInterfaceFluxes(isSubStep)) {
231 this->captureLocalFluxData();
232 }
233
234 if (this->collectOnIORank_.isParallel()) {
236
237 this->collectOnIORank_.collect({},
238 outputModule_->getBlockData(),
239 localWellData,
240 localWBP,
241 localGroupAndNetworkData,
242 localAquiferData,
243 localWellTestState,
244 this->outputModule_->getInterRegFlows(),
245 {},
246 {});
247
248 if (this->collectOnIORank_.isIORank()) {
249 auto& iregFlows = this->collectOnIORank_.globalInterRegFlows();
250
251 if (! iregFlows.readIsConsistent()) {
252 throw std::runtime_error {
253 "Inconsistent inter-region flow "
254 "region set names in parallel"
255 };
256 }
257
258 iregFlows.compress();
259 }
260
261 OPM_END_PARALLEL_TRY_CATCH("Collect to I/O rank: ",
262 this->simulator_.vanguard().grid().comm());
263 }
264
265
266 std::map<std::string, double> miscSummaryData;
267 std::map<std::string, std::vector<double>> regionData;
268 Inplace inplace;
269
270 {
271 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
272
273 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
274
275 if (this->collectOnIORank_.isIORank()){
276 inplace_ = inplace;
277 }
278 }
279
280 // Add TCPU
281 if (totalCpuTime != 0.0) {
282 miscSummaryData["TCPU"] = totalCpuTime;
283 }
285 miscSummaryData["NEWTON"] = this->sub_step_report_.total_newton_iterations;
286 }
288 miscSummaryData["MLINEARS"] = this->sub_step_report_.total_linear_iterations;
289 }
291 miscSummaryData["NLINEARS"] = static_cast<float>(this->sub_step_report_.total_linear_iterations) / this->sub_step_report_.total_newton_iterations;
292 }
293 if (this->sub_step_report_.min_linear_iterations != std::numeric_limits<unsigned int>::max()) {
294 miscSummaryData["NLINSMIN"] = this->sub_step_report_.min_linear_iterations;
295 }
297 miscSummaryData["NLINSMAX"] = this->sub_step_report_.max_linear_iterations;
298 }
300 miscSummaryData["MSUMLINS"] = this->simulation_report_.success.total_linear_iterations;
301 }
303 miscSummaryData["MSUMNEWT"] = this->simulation_report_.success.total_newton_iterations;
304 }
305
306 {
307 OPM_TIMEBLOCK(evalSummary);
308
309 const auto& blockData = this->collectOnIORank_.isParallel()
311 : this->outputModule_->getBlockData();
312
313 const auto& interRegFlows = this->collectOnIORank_.isParallel()
315 : this->outputModule_->getInterRegFlows();
316
317 this->evalSummary(reportStepNum,
318 curTime,
319 localWellData,
320 localWBP,
321 localGroupAndNetworkData,
322 localAquiferData,
323 blockData,
324 miscSummaryData,
325 regionData,
326 inplace,
327 this->outputModule_->initialInplace(),
328 interRegFlows,
329 this->summaryState(),
330 this->udqState());
331 }
332 }
333
336 {
337 const auto& fip = simulator_.vanguard().eclState().getEclipseConfig().fip();
338 if (!fip.output(FIPConfig::OutputField::FIELD) &&
339 !fip.output(FIPConfig::OutputField::RESV)) {
340 return;
341 }
342
343 const auto& gridView = simulator_.vanguard().gridView();
344 const int num_interior = detail::
346
347 this->outputModule_->
348 allocBuffers(num_interior, 0, false, false, /*isRestart*/ false);
349
350#ifdef _OPENMP
351#pragma omp parallel for
352#endif
353 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
354 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
355 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
356
357 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
358 }
359
360 std::map<std::string, double> miscSummaryData;
361 std::map<std::string, std::vector<double>> regionData;
362 Inplace inplace;
363 {
364 OPM_TIMEBLOCK(outputFipLogAndFipresvLog);
365
366 boost::posix_time::ptime start_time = boost::posix_time::from_time_t(simulator_.vanguard().schedule().getStartTime());
367
368 inplace = outputModule_->calc_inplace(miscSummaryData, regionData, simulator_.gridView().comm());
369
370 if (this->collectOnIORank_.isIORank()){
371 inplace_ = inplace;
372
373 outputModule_->outputFipAndResvLog(inplace_, 0, 0.0, start_time,
374 false, simulator_.gridView().comm());
375 }
376 }
377 }
378
379 void writeReports(const SimulatorTimer& timer) {
380 auto rstep = timer.reportStepNum();
381
382 if ((rstep > 0) && (this->collectOnIORank_.isIORank())){
383
384 const auto& rpt = this->schedule_[rstep-1].rpt_config.get();
385 if (rpt.contains("WELLS") && rpt.at("WELLS") > 0) {
386 outputModule_->outputTimeStamp("WELLS", timer.simulationTimeElapsed(), rstep, timer.currentDateTime());
387 outputModule_->outputProdLog(rstep-1);
388 outputModule_->outputInjLog(rstep-1);
389 outputModule_->outputCumLog(rstep-1);
390 }
391
392 outputModule_->outputFipAndResvLog(inplace_, rstep, timer.simulationTimeElapsed(),
393 timer.currentDateTime(), false, simulator_.gridView().comm());
394
395
396 OpmLog::note(""); // Blank line after all reports.
397 }
398 }
399
400 void writeOutput(data::Solution&& localCellData, bool isSubStep)
401 {
402 OPM_TIMEBLOCK(writeOutput);
403
404 const int reportStepNum = simulator_.episodeIndex() + 1;
405 this->prepareLocalCellData(isSubStep, reportStepNum);
406 this->outputModule_->outputErrorLog(simulator_.gridView().comm());
407
408 // output using eclWriter if enabled
409 auto localWellData = simulator_.problem().wellModel().wellData();
410 auto localGroupAndNetworkData = simulator_.problem().wellModel()
411 .groupAndNetworkData(reportStepNum);
412
413 auto localAquiferData = simulator_.problem().aquiferModel().aquiferData();
414 auto localWellTestState = simulator_.problem().wellModel().wellTestState();
415
416 const bool isFlowsn = this->outputModule_->hasFlowsn();
417 auto flowsn = this->outputModule_->getFlowsn();
418
419 const bool isFloresn = this->outputModule_->hasFloresn();
420 auto floresn = this->outputModule_->getFloresn();
421
422 // data::Solution localCellData = {};
423 if (! isSubStep || Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
424
425 if (localCellData.empty()) {
426 this->outputModule_->assignToSolution(localCellData);
427 }
428
429 // Add cell data to perforations for RFT output
430 this->outputModule_->addRftDataToWells(localWellData, reportStepNum);
431 }
432
433 if (this->collectOnIORank_.isParallel() ||
434 this->collectOnIORank_.doesNeedReordering())
435 {
436 // Note: We don't need WBP (well-block averaged pressures) or
437 // inter-region flow rate values in order to create restart file
438 // output. There's consequently no need to collect those
439 // properties on the I/O rank.
440
441 this->collectOnIORank_.collect(localCellData,
442 this->outputModule_->getBlockData(),
443 localWellData,
444 /* wbpData = */ {},
445 localGroupAndNetworkData,
446 localAquiferData,
447 localWellTestState,
448 /* interRegFlows = */ {},
449 flowsn,
450 floresn);
451 if (this->collectOnIORank_.isIORank()) {
452 this->outputModule_->assignGlobalFieldsToSolution(this->collectOnIORank_.globalCellData());
453 }
454 } else {
455 this->outputModule_->assignGlobalFieldsToSolution(localCellData);
456 }
457
458 if (this->collectOnIORank_.isIORank()) {
459 const Scalar curTime = simulator_.time() + simulator_.timeStepSize();
460 const Scalar nextStepSize = simulator_.problem().nextTimeStepSize();
461 std::optional<int> timeStepIdx;
462 if (Parameters::Get<Parameters::EnableWriteAllSolutions>()) {
463 timeStepIdx = simulator_.timeStepIndex();
464 }
465 this->doWriteOutput(reportStepNum, timeStepIdx, isSubStep,
466 std::move(localCellData),
467 std::move(localWellData),
468 std::move(localGroupAndNetworkData),
469 std::move(localAquiferData),
470 std::move(localWellTestState),
471 this->actionState(),
472 this->udqState(),
473 this->summaryState(),
474 this->simulator_.problem().thresholdPressure().getRestartVector(),
475 curTime, nextStepSize,
476 Parameters::Get<Parameters::EclOutputDoublePrecision>(),
477 isFlowsn, std::move(flowsn),
478 isFloresn, std::move(floresn));
479 }
480 }
481
483 {
484 bool enablePCHysteresis = simulator_.problem().materialLawManager()->enablePCHysteresis();
485 bool enableNonWettingHysteresis = simulator_.problem().materialLawManager()->enableNonWettingHysteresis();
486 bool enableWettingHysteresis = simulator_.problem().materialLawManager()->enableWettingHysteresis();
487 bool oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
488 bool gasActive = FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
489 bool waterActive = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
490 bool enableSwatinit = simulator_.vanguard().eclState().fieldProps().has_double("SWATINIT");
491 bool opm_rst_file = Parameters::Get<Parameters::EnableOpmRstFile>();
492 bool read_temp = enableEnergy || (opm_rst_file && enableTemperature);
493 std::vector<RestartKey> solutionKeys{
494 {"PRESSURE", UnitSystem::measure::pressure},
495 {"SWAT", UnitSystem::measure::identity, static_cast<bool>(waterActive)},
496 {"SGAS", UnitSystem::measure::identity, static_cast<bool>(gasActive)},
497 {"TEMP" , UnitSystem::measure::temperature, read_temp},
498 {"SSOLVENT" , UnitSystem::measure::identity, enableSolvent},
499 {"RS", UnitSystem::measure::gas_oil_ratio, FluidSystem::enableDissolvedGas()},
500 {"RV", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedOil()},
501 {"RVW", UnitSystem::measure::oil_gas_ratio, FluidSystem::enableVaporizedWater()},
502 {"SGMAX", UnitSystem::measure::identity, (enableNonWettingHysteresis && oilActive && gasActive)},
503 {"SHMAX", UnitSystem::measure::identity, (enableWettingHysteresis && oilActive && gasActive)},
504 {"SOMAX", UnitSystem::measure::identity, (enableNonWettingHysteresis && oilActive && waterActive) || simulator_.problem().vapparsActive(simulator_.episodeIndex())},
505 {"SOMIN", UnitSystem::measure::identity, (enablePCHysteresis && oilActive && gasActive)},
506 {"SWHY1", UnitSystem::measure::identity, (enablePCHysteresis && oilActive && waterActive)},
507 {"SWMAX", UnitSystem::measure::identity, (enableWettingHysteresis && oilActive && waterActive)},
508 {"PPCW", UnitSystem::measure::pressure, enableSwatinit}
509 };
510
511 const auto& inputThpres = eclState().getSimulationConfig().getThresholdPressure();
512 std::vector<RestartKey> extraKeys = {{"OPMEXTRA", UnitSystem::measure::identity, false},
513 {"THRESHPR", UnitSystem::measure::pressure, inputThpres.active()}};
514
515 {
516 const auto& tracers = simulator_.vanguard().eclState().tracer();
517 for (const auto& tracer : tracers) {
518 bool enableSolTracer = (tracer.phase == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
519 (tracer.phase == Phase::OIL && FluidSystem::enableVaporizedOil());
520 solutionKeys.emplace_back(tracer.fname(), UnitSystem::measure::identity, true);
521 solutionKeys.emplace_back(tracer.sname(), UnitSystem::measure::identity, enableSolTracer);
522 }
523 }
524
525 // The episodeIndex is rewined one back before beginRestart is called
526 // and can not be used here.
527 // We just ask the initconfig directly to be sure that we use the correct
528 // index.
529 const auto& initconfig = simulator_.vanguard().eclState().getInitConfig();
530 int restartStepIdx = initconfig.getRestartStep();
531
532 const auto& gridView = simulator_.vanguard().gridView();
533 unsigned numElements = gridView.size(/*codim=*/0);
534 outputModule_->allocBuffers(numElements, restartStepIdx, /*isSubStep=*/false, /*log=*/false, /*isRestart*/ true);
535
536 {
537 SummaryState& summaryState = simulator_.vanguard().summaryState();
538 Action::State& actionState = simulator_.vanguard().actionState();
539 auto restartValues = loadParallelRestart(this->eclIO_.get(), actionState, summaryState, solutionKeys, extraKeys,
540 gridView.grid().comm());
541 for (unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
542 unsigned globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
543 outputModule_->setRestart(restartValues.solution, elemIdx, globalIdx);
544 }
545
546 auto& tracer_model = simulator_.problem().tracerModel();
547 for (int tracer_index = 0; tracer_index < tracer_model.numTracers(); tracer_index++) {
548 // Free tracers
549 const auto& free_tracer_name = tracer_model.fname(tracer_index);
550 const auto& free_tracer_solution = restartValues.solution.template data<double>(free_tracer_name);
551 for (unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
552 unsigned globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
553 tracer_model.setFreeTracerConcentration(tracer_index, globalIdx, free_tracer_solution[globalIdx]);
554 }
555 // Solution tracer (only if DISGAS/VAPOIL are active for gas/oil tracers)
556 if ((tracer_model.phase(tracer_index) == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
557 (tracer_model.phase(tracer_index) == Phase::OIL && FluidSystem::enableVaporizedOil())) {
558 tracer_model.setEnableSolTracers(tracer_index, true);
559 const auto& sol_tracer_name = tracer_model.sname(tracer_index);
560 const auto& sol_tracer_solution = restartValues.solution.template data<double>(sol_tracer_name);
561 for (unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
562 unsigned globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
563 tracer_model.setSolTracerConcentration(tracer_index, globalIdx, sol_tracer_solution[globalIdx]);
564 }
565 }
566 else {
567 tracer_model.setEnableSolTracers(tracer_index, false);
568 for (unsigned elemIdx = 0; elemIdx < numElements; ++elemIdx) {
569 unsigned globalIdx = this->collectOnIORank_.localIdxToGlobalIdx(elemIdx);
570 tracer_model.setSolTracerConcentration(tracer_index, globalIdx, 0.0);
571 }
572 }
573 }
574
575 if (inputThpres.active()) {
576 Simulator& mutableSimulator = const_cast<Simulator&>(simulator_);
577 auto& thpres = mutableSimulator.problem().thresholdPressure();
578 const auto& thpresValues = restartValues.getExtra("THRESHPR");
579 thpres.setFromRestart(thpresValues);
580 }
581 restartTimeStepSize_ = restartValues.getExtra("OPMEXTRA")[0];
582 if (restartTimeStepSize_ <= 0)
583 restartTimeStepSize_ = std::numeric_limits<double>::max();
584
585 // initialize the well model from restart values
586 simulator_.problem().wellModel().initFromRestartFile(restartValues);
587
588 if (!restartValues.aquifer.empty())
589 simulator_.problem().mutableAquiferModel().initFromRestart(restartValues.aquifer);
590 }
591 }
592
594 {}
595
597 { return *outputModule_; }
598
600 { return *outputModule_; }
601
602 Scalar restartTimeStepSize() const
603 { return restartTimeStepSize_; }
604
605 template <class Serializer>
606 void serializeOp(Serializer& serializer)
607 {
608 serializer(*outputModule_);
609 }
610
611private:
612 static bool enableEclOutput_()
613 {
614 static bool enable = Parameters::Get<Parameters::EnableEclOutput>();
615 return enable;
616 }
617
618 const EclipseState& eclState() const
619 { return simulator_.vanguard().eclState(); }
620
621 SummaryState& summaryState()
622 { return simulator_.vanguard().summaryState(); }
623
624 Action::State& actionState()
625 { return simulator_.vanguard().actionState(); }
626
627 UDQState& udqState()
628 { return simulator_.vanguard().udqState(); }
629
630 const Schedule& schedule() const
631 { return simulator_.vanguard().schedule(); }
632
633 void prepareLocalCellData(const bool isSubStep,
634 const int reportStepNum)
635 {
636 OPM_TIMEBLOCK(prepareLocalCellData);
637
638 if (this->outputModule_->localDataValid()) {
639 return;
640 }
641
642 const auto& gridView = simulator_.vanguard().gridView();
643 const bool log = this->collectOnIORank_.isIORank();
644
645 const int num_interior = detail::
647 this->outputModule_->
648 allocBuffers(num_interior, reportStepNum,
649 isSubStep && !Parameters::Get<Parameters::EnableWriteAllSolutions>(),
650 log, /*isRestart*/ false);
651
652 ElementContext elemCtx(simulator_);
653
655
656 {
657 OPM_TIMEBLOCK(prepareCellBasedData);
658
659 this->outputModule_->prepareDensityAccumulation();
660
661 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
662 elemCtx.updatePrimaryStencil(elem);
663 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
664
665 this->outputModule_->processElement(elemCtx);
666 }
667
668 this->outputModule_->accumulateDensityParallel();
669 }
670
671 if constexpr (enableMech) {
672 if (simulator_.vanguard().eclState().runspec().mech()) {
673 OPM_TIMEBLOCK(prepareMechData);
674 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
675 elemCtx.updatePrimaryStencil(elem);
676 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
677 outputModule_->processElementMech(elemCtx);
678 }
679 }
680 }
681
682 if (! this->simulator_.model().linearizer().getFlowsInfo().empty()) {
683 OPM_TIMEBLOCK(prepareFlowsData);
684 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
685 elemCtx.updatePrimaryStencil(elem);
686 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
687
688 this->outputModule_->processElementFlows(elemCtx);
689 }
690 }
691
692 {
693 OPM_TIMEBLOCK(prepareBlockData);
694 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
695 elemCtx.updatePrimaryStencil(elem);
696 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
697
698 this->outputModule_->processElementBlockData(elemCtx);
699 }
700 }
701
702 {
703 OPM_TIMEBLOCK(prepareFluidInPlace);
704
705#ifdef _OPENMP
706#pragma omp parallel for
707#endif
708 for (int dofIdx = 0; dofIdx < num_interior; ++dofIdx) {
709 const auto& intQuants = *simulator_.model().cachedIntensiveQuantities(dofIdx, /*timeIdx=*/0);
710 const auto totVolume = simulator_.model().dofTotalVolume(dofIdx);
711
712 this->outputModule_->updateFluidInPlace(dofIdx, intQuants, totVolume);
713 }
714 }
715
716 this->outputModule_->validateLocalData();
717
718 OPM_END_PARALLEL_TRY_CATCH("EclWriter::prepareLocalCellData() failed: ",
719 this->simulator_.vanguard().grid().comm());
720 }
721
722 void captureLocalFluxData()
723 {
724 OPM_TIMEBLOCK(captureLocalData);
725
726 const auto& gridView = this->simulator_.vanguard().gridView();
727 const auto timeIdx = 0u;
728
729 auto elemCtx = ElementContext { this->simulator_ };
730
731 const auto elemMapper = ElementMapper { gridView, Dune::mcmgElementLayout() };
732 const auto activeIndex = [&elemMapper](const Element& e)
733 {
734 return elemMapper.index(e);
735 };
736
737 const auto cartesianIndex = [this](const int elemIndex)
738 {
739 return this->cartMapper_.cartesianIndex(elemIndex);
740 };
741
742 this->outputModule_->initializeFluxData();
743
745
746 for (const auto& elem : elements(gridView, Dune::Partitions::interiorBorder)) {
747 elemCtx.updateStencil(elem);
748 elemCtx.updateIntensiveQuantities(timeIdx);
749 elemCtx.updateExtensiveQuantities(timeIdx);
750
751 this->outputModule_->processFluxes(elemCtx, activeIndex, cartesianIndex);
752 }
753
754 OPM_END_PARALLEL_TRY_CATCH("EclWriter::captureLocalFluxData() failed: ",
755 this->simulator_.vanguard().grid().comm())
756
757 this->outputModule_->finalizeFluxData();
758 }
759
760 Simulator& simulator_;
761 std::unique_ptr<OutputBlackOilModule<TypeTag> > outputModule_;
762 Scalar restartTimeStepSize_;
763 int rank_ ;
764 Inplace inplace_;
765};
766
767} // namespace Opm
768
769#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
void collect(const data::Solution &localCellData, const std::map< std::pair< std::string, int >, double > &localBlockData, const data::Wells &localWellData, const data::WellBlockAveragePressures &localWBPData, const data::GroupAndNetworkValues &localGroupAndNetworkData, const data::Aquifers &localAquiferData, const WellTestState &localWellTestState, const InterRegFlowMap &interRegFlows, const std::array< FlowsData< double >, 3 > &localFlowsn, const std::array< FlowsData< double >, 3 > &localFloresn)
Definition: CollectDataOnIORank_impl.hpp:969
int localIdxToGlobalIdx(unsigned localIdx) const
Definition: CollectDataOnIORank_impl.hpp:1091
InterRegFlowMap & globalInterRegFlows()
Definition: CollectDataOnIORank.hpp:113
bool isParallel() const
Definition: CollectDataOnIORank.hpp:128
bool isIORank() const
Definition: CollectDataOnIORank.hpp:125
const std::map< std::pair< std::string, int >, double > & globalBlockData() const
Definition: CollectDataOnIORank.hpp:89
const data::Solution & globalCellData() const
Definition: CollectDataOnIORank.hpp:92
Definition: EclGenericWriter.hpp:65
void doWriteOutput(const int reportStepNum, const std::optional< int > timeStepNum, const bool isSubStep, data::Solution &&localCellData, data::Wells &&localWellData, data::GroupAndNetworkValues &&localGroupAndNetworkData, data::Aquifers &&localAquiferData, WellTestState &&localWTestState, const Action::State &actionState, const UDQState &udqState, const SummaryState &summaryState, const std::vector< GetPropType< TypeTag, Properties::Scalar > > &thresholdPressure, GetPropType< TypeTag, Properties::Scalar > curTime, GetPropType< TypeTag, Properties::Scalar > nextStepSize, bool doublePrecision, bool isFlowsn, std::array< FlowsData< double >, 3 > &&flowsn, bool isFloresn, std::array< FlowsData< double >, 3 > &&floresn)
Definition: EclGenericWriter_impl.hpp:531
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:624
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:101
void writeInitialFIPReport()
Writes the initial FIP report as configured in RPTSOL.
Definition: EclWriter.hpp:335
void endRestart()
Definition: EclWriter.hpp:593
const OutputBlackOilModule< TypeTag > & outputModule() const
Definition: EclWriter.hpp:596
Scalar restartTimeStepSize() const
Definition: EclWriter.hpp:602
~EclWriter()
Definition: EclWriter.hpp:178
void writeOutput(data::Solution &&localCellData, bool isSubStep)
Definition: EclWriter.hpp:400
EclWriter(Simulator &simulator)
Definition: EclWriter.hpp:138
static void registerParameters()
Definition: EclWriter.hpp:124
void writeReports(const SimulatorTimer &timer)
Definition: EclWriter.hpp:379
const EquilGrid & globalGrid() const
Definition: EclWriter.hpp:181
void beginRestart()
Definition: EclWriter.hpp:482
OutputBlackOilModule< TypeTag > & mutableOutputModule() const
Definition: EclWriter.hpp:599
void evalSummaryState(bool isSubStep)
collect and pass data and pass it to eclIO writer
Definition: EclWriter.hpp:189
void serializeOp(Serializer &serializer)
Definition: EclWriter.hpp:606
Output module for the results black oil model writing in ECL binary format.
Definition: OutputBlackoilModule.hpp:90
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: OutputBlackoilModule.hpp:182
virtual int reportStepNum() const
Current report step number. This might differ from currentStepNum in case of sub stepping.
Definition: SimulatorTimerInterface.hpp:50
Definition: SimulatorTimer.hpp:39
virtual boost::posix_time::ptime currentDateTime() const
Return the current time as a posix time object.
double simulationTimeElapsed() const override
Definition: blackoilnewtonmethodparams.hpp:31
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition: parametersystem.hpp:185
std::size_t countLocalInteriorCellsGridView(const GridView &gridView)
Get the number of local interior cells in a grid view.
Definition: countGlobalCells.hpp:47
Definition: blackoilboundaryratevector.hh:37
void eclBroadcast(Parallel::Communication, T &)
RestartValue loadParallelRestart(const EclipseIO *eclIO, Action::State &actionState, SummaryState &summaryState, const std::vector< RestartKey > &solutionKeys, const std::vector< RestartKey > &extraKeys, Parallel::Communication comm)
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
Definition: EclWriter.hpp:63
static constexpr bool value
Definition: EclWriter.hpp:63
Definition: EclWriter.hpp:60
static constexpr bool value
Definition: EclWriter.hpp:60
Definition: EclWriter.hpp:57
static constexpr bool value
Definition: EclWriter.hpp:57
Definition: EclWriter.hpp:70
static constexpr bool value
Definition: EclWriter.hpp:70
Definition: EclWriter.hpp:67
static constexpr bool value
Definition: EclWriter.hpp:67
SimulatorReportSingle success
Definition: SimulatorReport.hpp:101
unsigned int min_linear_iterations
Definition: SimulatorReport.hpp:51
unsigned int total_newton_iterations
Definition: SimulatorReport.hpp:49
unsigned int max_linear_iterations
Definition: SimulatorReport.hpp:52
unsigned int total_linear_iterations
Definition: SimulatorReport.hpp:50