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