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