FlowProblemBlackoil.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 Copyright 2023 INRIA
5 Copyright 2024 SINTEF Digital
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 2 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21
22 Consult the COPYING file in the top-level source directory of this
23 module for the precise wording of the license and the list of
24 copyright holders.
25*/
31#ifndef OPM_FLOW_PROBLEM_BLACK_HPP
32#define OPM_FLOW_PROBLEM_BLACK_HPP
33
34#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
35#include <opm/material/fluidsystems/blackoilpvt/DryGasPvt.hpp>
36#include <opm/material/fluidsystems/blackoilpvt/WetGasPvt.hpp>
37#include <opm/material/fluidsystems/blackoilpvt/LiveOilPvt.hpp>
38#include <opm/material/fluidsystems/blackoilpvt/DeadOilPvt.hpp>
39#include <opm/material/fluidsystems/blackoilpvt/ConstantCompressibilityOilPvt.hpp>
40#include <opm/material/fluidsystems/blackoilpvt/ConstantCompressibilityWaterPvt.hpp>
41
44
45#include <opm/output/eclipse/EclipseIO.hpp>
46
47#include <opm/input/eclipse/Units/Units.hpp>
48
58
60
61#if HAVE_DAMARIS
63#endif
64
65#include <algorithm>
66#include <cstddef>
67#include <functional>
68#include <limits>
69#include <memory>
70#include <stdexcept>
71#include <string>
72#include <string_view>
73#include <vector>
74
75namespace Opm {
76
83template <class TypeTag>
84class FlowProblemBlackoil : public FlowProblem<TypeTag>
85{
86 // TODO: the naming of the Types might be able to be adjusted
87public:
89
90private:
91 using typename FlowProblemType::Scalar;
92 using typename FlowProblemType::Simulator;
93 using typename FlowProblemType::GridView;
94 using typename FlowProblemType::FluidSystem;
95 using typename FlowProblemType::Vanguard;
97 using typename FlowProblemType::EqVector;
103
104 // TODO: potentially some cleaning up depending on the usage later here
120
124
128
130 using typename FlowProblemType::RateVector;
132 using typename FlowProblemType::Indices;
134 using typename FlowProblemType::ElementContext;
135
136 using typename FlowProblemType::MaterialLaw;
137 using typename FlowProblemType::DimMatrix;
138
139 enum { enableDissolvedGas = Indices::compositionSwitchIdx >= 0 };
140 enum { enableVapwat = getPropValue<TypeTag, Properties::EnableVapwat>() };
141 enum { enableDisgasInWater = getPropValue<TypeTag, Properties::EnableDisgasInWater>() };
142
143 using SolventModule = BlackOilSolventModule<TypeTag>;
144 using PolymerModule = BlackOilPolymerModule<TypeTag>;
145 using FoamModule = BlackOilFoamModule<TypeTag>;
146 using BrineModule = BlackOilBrineModule<TypeTag>;
147 using ExtboModule = BlackOilExtboModule<TypeTag>;
148 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag>;
149 using DispersionModule = BlackOilDispersionModule<TypeTag, enableDispersion>;
150 using DiffusionModule = BlackOilDiffusionModule<TypeTag, enableDiffusion>;
151 using ConvectiveMixingModule = BlackOilConvectiveMixingModule<TypeTag, enableConvectiveMixing>;
152 using ModuleParams = BlackoilModuleParams<ConvectiveMixingModuleParam<Scalar>>;
153 using HybridNewton = BlackOilHybridNewton<TypeTag>;
154
155 using InitialFluidState = typename EquilInitializer<TypeTag>::ScalarFluidState;
156 using EclWriterType = EclWriter<TypeTag, OutputBlackOilModule<TypeTag> >;
157 using IndexTraits = typename FluidSystem::IndexTraitsType;
158#if HAVE_DAMARIS
159 using DamarisWriterType = DamarisWriter<TypeTag>;
160#endif
161
162
163public:
166
170 static void registerParameters()
171 {
173
175#if HAVE_DAMARIS
176 DamarisWriterType::registerParameters();
177#endif
179 }
180
184 explicit FlowProblemBlackoil(Simulator& simulator)
185 : FlowProblemType(simulator)
186 , thresholdPressures_(simulator)
187 , mixControls_(simulator.vanguard().schedule())
188 , actionHandler_(simulator.vanguard().eclState(),
189 simulator.vanguard().schedule(),
190 simulator.vanguard().actionState(),
191 simulator.vanguard().summaryState(),
192 this->wellModel_,
193 simulator.vanguard().grid().comm())
194 , hybridNewton_(simulator)
195 {
196 this->model().addOutputModule(std::make_unique<VtkTracerModule<TypeTag>>(simulator));
197
198 // Tell the black-oil extensions to initialize their internal data structures
199 const auto& vanguard = simulator.vanguard();
200
201 BlackOilBrineParams<Scalar> brineParams;
202 brineParams.template initFromState<enableBrine,
203 enableSaltPrecipitation>(vanguard.eclState());
204 BrineModule::setParams(std::move(brineParams));
205
206 DiffusionModule::initFromState(vanguard.eclState());
207 DispersionModule::initFromState(vanguard.eclState());
208
209 BlackOilExtboParams<Scalar> extboParams;
210 extboParams.template initFromState<enableExtbo>(vanguard.eclState());
211 ExtboModule::setParams(std::move(extboParams));
212
214 foamParams.template initFromState<enableFoam>(vanguard.eclState());
215 FoamModule::setParams(std::move(foamParams));
216
217 BlackOilBioeffectsParams<Scalar> bioeffectsParams;
218 bioeffectsParams.template initFromState<enableBioeffects, enableMICP>(vanguard.eclState());
219 BioeffectsModule::setParams(std::move(bioeffectsParams));
220
221 BlackOilPolymerParams<Scalar> polymerParams;
222 polymerParams.template initFromState<enablePolymer, enablePolymerMolarWeight>(vanguard.eclState());
223 PolymerModule::setParams(std::move(polymerParams));
224
225 BlackOilSolventParams<Scalar> solventParams;
226 solventParams.template initFromState<enableSolvent>(vanguard.eclState(), vanguard.schedule());
227 SolventModule::setParams(std::move(solventParams));
228
229 // create the ECL writer
230 eclWriter_ = std::make_unique<EclWriterType>(simulator);
231 enableEclOutput_ = Parameters::Get<Parameters::EnableEclOutput>();
232
233#if HAVE_DAMARIS
234 // create Damaris writer
235 damarisWriter_ = std::make_unique<DamarisWriterType>(simulator);
236 enableDamarisOutput_ = Parameters::Get<Parameters::EnableDamarisOutput>();
237#endif
238 }
239
243 void beginEpisode() override
244 {
246
247 auto& simulator = this->simulator();
248
249 const int episodeIdx = simulator.episodeIndex();
250 const auto& schedule = simulator.vanguard().schedule();
251
252 // Evaluate UDQ assign statements to make sure the settings are
253 // available as UDA controls for the current report step.
254 this->actionHandler_
255 .evalUDQAssignments(episodeIdx, simulator.vanguard().udqState());
256
257 if (episodeIdx >= 0) {
258 const auto& oilVap = schedule[episodeIdx].oilvap();
259 if (oilVap.getType() == OilVaporizationProperties::OilVaporization::VAPPARS) {
260 FluidSystem::setVapPars(oilVap.vap1(), oilVap.vap2());
261 }
262 else {
263 FluidSystem::setVapPars(0.0, 0.0);
264 }
265 }
266
267 ConvectiveMixingModule::beginEpisode(simulator.vanguard().eclState(), schedule, episodeIdx,
268 this->moduleParams_.convectiveMixingModuleParam);
269 }
270
274 void beginTimeStep() override
275 {
278 }
279
284 {
285 // TODO: there should be room to remove duplication for this
286 // function, but there is relatively complicated logic in the
287 // function calls here. Some refactoring is needed.
288 FlowProblemType::finishInit();
289
290 auto& simulator = this->simulator();
291
292 auto finishTransmissibilities = [updated = false, this]() mutable
293 {
294 if (updated) { return; }
295
296 this->transmissibilities_.finishInit([&vg = this->simulator().vanguard()](const unsigned int it) {
297 return vg.gridIdxToEquilGridIdx(it);
298 });
299
300 updated = true;
301 };
302
303 // calculating the TRANX, TRANY, TRANZ and NNC for output purpose
304 // for parallel running, it is based on global trans_
305 // for serial running, it is based on the transmissibilities_
306 // we try to avoid for the parallel running, has both global trans_ and transmissibilities_ allocated at the same time
307 if (enableEclOutput_) {
308 if (simulator.vanguard().grid().comm().size() > 1) {
309 if (simulator.vanguard().grid().comm().rank() == 0)
310 eclWriter_->setTransmissibilities(&simulator.vanguard().globalTransmissibility());
311 } else {
312 finishTransmissibilities();
313 eclWriter_->setTransmissibilities(&simulator.problem().eclTransmissibilities());
314 }
315
316 std::function<unsigned int(unsigned int)> equilGridToGrid = [&simulator](unsigned int i) {
317 return simulator.vanguard().gridEquilIdxToGridIdx(i);
318 };
319
320 this->eclWriter_->extractOutputTransAndNNC(equilGridToGrid);
321 }
322 simulator.vanguard().releaseGlobalTransmissibilities();
323
324 const auto& eclState = simulator.vanguard().eclState();
325 const auto& schedule = simulator.vanguard().schedule();
326
327 // Set the start time of the simulation
328 simulator.setStartTime(schedule.getStartTime());
329 simulator.setEndTime(schedule.simTime(schedule.size() - 1));
330
331 // We want the episode index to be the same as the report step index to make
332 // things simpler, so we have to set the episode index to -1 because it is
333 // incremented by endEpisode(). The size of the initial time step and
334 // length of the initial episode is set to zero for the same reason.
335 simulator.setEpisodeIndex(-1);
336 simulator.setEpisodeLength(0.0);
337
338 // the "NOGRAV" keyword from Frontsim or setting the EnableGravity to false
339 // disables gravity, else the standard value of the gravity constant at sea level
340 // on earth is used
341 this->gravity_ = 0.0;
342 if (Parameters::Get<Parameters::EnableGravity>() &&
343 eclState.getInitConfig().hasGravity())
344 {
345 // unit::gravity is 9.80665 m^2/s--i.e., standard measure at Tellus equator.
346 this->gravity_[dim - 1] = unit::gravity;
347 }
348
349 if (this->enableTuning_) {
350 // if support for the TUNING keyword is enabled, we get the initial time
351 // steping parameters from it instead of from command line parameters
352 const auto& tuning = schedule[0].tuning();
353 this->initialTimeStepSize_ = tuning.TSINIT.has_value() ? tuning.TSINIT.value() : -1.0;
354 this->maxTimeStepAfterWellEvent_ = tuning.TMAXWC;
355 }
356
357 // conserve inner energy instead of enthalpy if TEMP is used
358 // or THERMAL and parameter ConserveInnerEnergyThermal is true (default false)
359 bool isThermal = eclState.getSimulationConfig().isThermal();
360 bool isTemp = eclState.getSimulationConfig().isTemp();
361 bool conserveInnerEnergy = isTemp || (isThermal && Parameters::Get<Parameters::ConserveInnerEnergyThermal>());
362 FluidSystem::setEnergyEqualEnthalpy(conserveInnerEnergy);
363
364 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
365 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
366 this->maxOilSaturation_.resize(this->model().numGridDof(), 0.0);
367 }
368
369 this->readRockParameters_(simulator.vanguard().cellCenterDepths(),
370 [&simulator](const unsigned idx)
371 {
372 std::array<int,dim> coords;
373 simulator.vanguard().cartesianCoordinate(idx, coords);
374 std::ranges::transform(coords, coords.begin(),
375 [](const auto c) { return c + 1; });
376 return coords;
377 });
378
381
382 // write the static output files (EGRID, INIT)
383 if (enableEclOutput_) {
384 this->eclWriter_->writeInit();
385 }
386
387 finishTransmissibilities();
388
389 const auto& initconfig = eclState.getInitConfig();
390 this->tracerModel_.init(initconfig.restartRequested());
391 if (initconfig.restartRequested()) {
393 }
394 else {
395 this->readInitialCondition_();
396 }
397 this->temperatureModel_.init();
398 this->tracerModel_.prepareTracerBatches();
399
400 this->updatePffDofData_();
401
402 if constexpr (getPropValue<TypeTag, Properties::EnablePolymer>()) {
403 const auto& vanguard = this->simulator().vanguard();
404 const auto& gridView = vanguard.gridView();
405 const int numElements = gridView.size(/*codim=*/0);
406 this->polymer_.maxAdsorption.resize(numElements, 0.0);
407 }
408
410
411 // compute and set eq weights based on initial b values
413
415 this->drift_.resize(this->model().numGridDof());
416 this->drift_ = 0.0;
417 }
418
419 // after finishing the initialization and writing the initial solution, we move
420 // to the first "real" episode/report step
421 // for restart the episode index and start is already set
422 if (!initconfig.restartRequested() && !eclState.getIOConfig().initOnly()) {
423 simulator.startNextEpisode(schedule.seconds(1));
424 simulator.setEpisodeIndex(0);
425 simulator.setTimeStepIndex(0);
426 }
427
428 if (Parameters::Get<Parameters::CheckSatfuncConsistency>() &&
430 {
431 // User requested that saturation functions be checked for
432 // consistency and essential/critical requirements are not met.
433 // Abort simulation run.
434 //
435 // Note: We need synchronisation here lest ranks other than the
436 // I/O rank throw exceptions too early thereby risking an
437 // incomplete failure report being shown to the user.
438 this->simulator().vanguard().grid().comm().barrier();
439
440 throw std::domain_error {
441 "Saturation function end-points do not "
442 "meet requisite consistency conditions"
443 };
444 }
445
446 // TODO: move to the end for later refactoring of the function finishInit()
447 //
448 // deal with DRSDT
449 this->mixControls_.init(this->model().numGridDof(),
450 this->episodeIndex(),
451 eclState.runspec().tabdims().getNumPVTTables());
452
453 if (this->enableVtkOutput_() && eclState.getIOConfig().initOnly()) {
454 simulator.setTimeStepSize(0.0);
455 simulator.model().applyInitialSolution();
457 }
458
459 if (!eclState.getIOConfig().initOnly()) {
460 if (!this->enableTuning_ && eclState.getSimulationConfig().anyTUNING()) {
461 OpmLog::info("\nThe deck has TUNING in the SCHEDULE section, but "
462 "it is ignored due\nto the flag --enable-tuning=false. "
463 "Set this flag to true to activate it.\n"
464 "Manually tuning the simulator with the TUNING keyword may "
465 "increase run time.\nIt is recommended using the simulator's "
466 "default tuning (--enable-tuning=false).");
467 }
468 }
469 }
470
474 void endTimeStep() override
475 {
476 FlowProblemType::endTimeStep();
477 this->endStepApplyAction();
478 }
479
481 {
482 // After the solution is updated, the values in output module needs
483 // also updated.
484 this->eclWriter().mutableOutputModule().invalidateLocalData();
485
486 // For CpGrid with LGRs, ecl/vtk output is not supported yet.
487 const auto& grid = this->simulator().vanguard().gridView().grid();
488
489 using GridType = std::remove_cv_t<std::remove_reference_t<decltype(grid)>>;
490 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
491 if (!isCpGrid || (grid.maxLevel() == 0)) {
492 this->eclWriter_->evalSummaryState(!this->episodeWillBeOver());
493 }
494
495 {
496 OPM_TIMEBLOCK(applyActions);
497
498 const int episodeIdx = this->episodeIndex();
499 auto& simulator = this->simulator();
500
501 // Clear out any existing events as these have already been
502 // processed when we're running an action block
503 this->simulator().vanguard().schedule().clearEvents(episodeIdx);
504
505 // Re-ordering in case of Alugrid
506 this->actionHandler_
507 .applyActions(episodeIdx, simulator.time() + simulator.timeStepSize(),
508 [this](const bool global)
509 {
510 using TransUpdateQuantities = typename
511 Vanguard::TransmissibilityType::TransUpdateQuantities;
512
513 this->transmissibilities_
514 .update(global, TransUpdateQuantities::All,
515 [&vg = this->simulator().vanguard()]
516 (const unsigned int i)
517 {
518 return vg.gridIdxToEquilGridIdx(i);
519 });
520 });
521 }
522 }
523
527 void endEpisode() override
528 {
529 OPM_TIMEBLOCK(endEpisode);
530
531 // Rerun UDQ assignents following action processing on the final
532 // time step of this episode to make sure that any UDQ ASSIGN
533 // operations triggered in action blocks take effect. This is
534 // mainly to work around a shortcoming of the ScheduleState copy
535 // constructor which clears pending UDQ assignments under the
536 // assumption that all such assignments have been processed. If an
537 // action block happens to trigger on the final time step of an
538 // episode and that action block runs a UDQ assignment, then that
539 // assignment would be dropped and the rest of the simulator will
540 // never see its effect without this hack.
541 this->actionHandler_
542 .evalUDQAssignments(this->episodeIndex(), this->simulator().vanguard().udqState());
543
544 FlowProblemType::endEpisode();
545 }
546
547 void writeReports(const SimulatorTimer& timer)
548 {
549 if (this->enableEclOutput_) {
550 this->eclWriter_->writeReports(timer);
551 }
552 }
553
554
559 void writeOutput(const bool verbose) override
560 {
561 FlowProblemType::writeOutput(verbose);
562
563 const auto isSubStep = !this->episodeWillBeOver();
564
565 auto localCellData = data::Solution {};
566
567#if HAVE_DAMARIS
568 // N.B. the Damaris output has to be done before the ECL output as the ECL one
569 // does all kinds of std::move() relocation of data
570 if (this->enableDamarisOutput_ && (this->damarisWriter_ != nullptr)) {
571 this->damarisWriter_->writeOutput(localCellData, isSubStep);
572 }
573#endif
574
575 if (this->enableEclOutput_ && (this->eclWriter_ != nullptr)) {
576 this->eclWriter_->writeOutput(std::move(localCellData), isSubStep);
577 }
578 }
579
581 {
582 OPM_TIMEBLOCK(finalizeOutput);
583 // this will write all pending output to disk
584 // to avoid corruption of output files
585 eclWriter_.reset();
586 }
587
588
593 {
594 FlowProblemType::initialSolutionApplied();
595
596 // let the object for threshold pressures initialize itself. this is done only at
597 // this point, because determining the threshold pressures may require to access
598 // the initial solution.
599 this->thresholdPressures_.finishInit();
600
601 // For CpGrid with LGRs, ecl-output is not supported yet.
602 const auto& grid = this->simulator().vanguard().gridView().grid();
603
604 using GridType = std::remove_cv_t<std::remove_reference_t<decltype(grid)>>;
605 constexpr bool isCpGrid = std::is_same_v<GridType, Dune::CpGrid>;
606 // Skip - for now - calculate the initial fip values for CpGrid with LGRs.
607 if (!isCpGrid || (grid.maxLevel() == 0)) {
608 if (this->simulator().episodeIndex() == 0) {
609 eclWriter_->writeInitialFIPReport();
610 }
611 }
612 }
613
615 unsigned globalDofIdx,
616 unsigned timeIdx) const override
617 {
618 this->aquiferModel_.addToSource(rate, globalDofIdx, timeIdx);
619
620 // Add source term from deck
621 const auto& source = this->simulator().vanguard().schedule()[this->episodeIndex()].source();
622 std::array<int,3> ijk;
623 this->simulator().vanguard().cartesianCoordinate(globalDofIdx, ijk);
624
625 if (source.hasSource(ijk)) {
626 const int pvtRegionIdx = this->pvtRegionIndex(globalDofIdx);
627 static std::array<SourceComponent, 3> sc_map = {SourceComponent::WATER, SourceComponent::OIL, SourceComponent::GAS};
628 static std::array<int, 3> phidx_map = {FluidSystem::waterPhaseIdx, FluidSystem::oilPhaseIdx, FluidSystem::gasPhaseIdx};
629 static std::array<int, 3> cidx_map = {waterCompIdx, oilCompIdx, gasCompIdx};
630
631 for (unsigned i = 0; i < phidx_map.size(); ++i) {
632 const auto phaseIdx = phidx_map[i];
633 const auto sourceComp = sc_map[i];
634 const auto compIdx = cidx_map[i];
635 if (!FluidSystem::phaseIsActive(phaseIdx)) {
636 continue;
637 }
638 Scalar mass_rate = source.rate(ijk, sourceComp) / this->model().dofTotalVolume(globalDofIdx);
639 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
640 mass_rate /= FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
641 }
642 rate[FluidSystem::canonicalToActiveCompIdx(compIdx)] += mass_rate;
643 }
644
645 if constexpr (enableSolvent) {
646 Scalar mass_rate = source.rate(ijk, SourceComponent::SOLVENT) / this->model().dofTotalVolume(globalDofIdx);
647 if constexpr (getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>()) {
648 const auto& solventPvt = SolventModule::solventPvt();
649 mass_rate /= solventPvt.referenceDensity(pvtRegionIdx);
650 }
651 rate[Indices::contiSolventEqIdx] += mass_rate;
652 }
653 if constexpr (enablePolymer) {
654 rate[Indices::polymerConcentrationIdx] += source.rate(ijk, SourceComponent::POLYMER) / this->model().dofTotalVolume(globalDofIdx);
655 }
656 if constexpr (enableMICP) {
657 rate[Indices::microbialConcentrationIdx] += source.rate(ijk, SourceComponent::MICR) / this->model().dofTotalVolume(globalDofIdx);
658 rate[Indices::oxygenConcentrationIdx] += source.rate(ijk, SourceComponent::OXYG) / this->model().dofTotalVolume(globalDofIdx);
659 rate[Indices::ureaConcentrationIdx] += source.rate(ijk, SourceComponent::UREA) / (this->model().dofTotalVolume(globalDofIdx));
660 }
661 if constexpr (energyModuleType == EnergyModules::FullyImplicitThermal) {
662 for (unsigned i = 0; i < phidx_map.size(); ++i) {
663 const auto phaseIdx = phidx_map[i];
664 if (!FluidSystem::phaseIsActive(phaseIdx)) {
665 continue;
666 }
667 const auto sourceComp = sc_map[i];
668 const auto source_hrate = source.hrate(ijk, sourceComp);
669 if (source_hrate) {
670 rate[Indices::contiEnergyEqIdx] += source_hrate.value() / this->model().dofTotalVolume(globalDofIdx);
671 } else {
672 const auto& intQuants = this->simulator().model().intensiveQuantities(globalDofIdx, /*timeIdx*/ 0);
673 auto fs = intQuants.fluidState();
674 // if temperature is not set, use cell temperature as default
675 const auto source_temp = source.temperature(ijk, sourceComp);
676 if (source_temp) {
677 Scalar temperature = source_temp.value();
678 fs.setTemperature(temperature);
679 }
680 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, pvtRegionIdx);
681 Scalar mass_rate = source.rate(ijk, sourceComp)/ this->model().dofTotalVolume(globalDofIdx);
682 Scalar energy_rate = getValue(h)*mass_rate;
683 rate[Indices::contiEnergyEqIdx] += energy_rate;
684 }
685 }
686 }
687 }
688
689 // if requested, compensate systematic mass loss for cells which were "well
690 // behaved" in the last time step
691 if (this->enableDriftCompensation_) {
692 const auto& simulator = this->simulator();
693 const auto& model = this->model();
694
695 // we use a lower tolerance for the compensation too
696 // assure the added drift from the last step does not
697 // cause convergence issues on the current step
698 Scalar maxCompensation = model.newtonMethod().tolerance()/10;
699 Scalar poro = this->porosity(globalDofIdx, timeIdx);
700 Scalar dt = simulator.timeStepSize();
701 EqVector dofDriftRate = this->drift_[globalDofIdx];
702 dofDriftRate /= dt*model.dofTotalVolume(globalDofIdx);
703
704 // restrict drift compensation to the CNV tolerance
705 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
706 Scalar cnv = std::abs(dofDriftRate[eqIdx])*dt*model.eqWeight(globalDofIdx, eqIdx)/poro;
707 if (cnv > maxCompensation) {
708 dofDriftRate[eqIdx] *= maxCompensation/cnv;
709 }
710 }
711
712 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx)
713 rate[eqIdx] -= dofDriftRate[eqIdx];
714 }
715 }
716
720 template <class LhsEval, class Callback>
721 LhsEval permFactTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx, Callback& obtain) const
722 {
723 OPM_TIMEBLOCK_LOCAL(permFactTransMultiplier, Subsystem::PvtProps);
724 if constexpr (enableSaltPrecipitation) {
725 const auto& fs = intQuants.fluidState();
726 unsigned tableIdx = this->simulator().problem().satnumRegionIndex(elementIdx);
727 LhsEval porosityFactor = obtain(1. - fs.saltSaturation());
728 porosityFactor = min(porosityFactor, 1.0);
729 const auto& permfactTable = BrineModule::permfactTable(tableIdx);
730 return permfactTable.eval(porosityFactor, /*extrapolation=*/true);
731 }
732 else if constexpr (enableBioeffects) {
733 return obtain(intQuants.permFactor());
734 }
735 else {
736 return 1.0;
737 }
738 }
739
740 // temporary solution to facilitate output of initial state from flow
741 const InitialFluidState& initialFluidState(unsigned globalDofIdx) const
742 { return initialFluidStates_[globalDofIdx]; }
743
744 std::vector<InitialFluidState>& initialFluidStates()
745 { return initialFluidStates_; }
746
747 const std::vector<InitialFluidState>& initialFluidStates() const
748 { return initialFluidStates_; }
749
750 const EclipseIO& eclIO() const
751 { return eclWriter_->eclIO(); }
752
754 { return eclWriter_->setSubStepReport(report); }
755
757 { return eclWriter_->setSimulationReport(report); }
758
759 InitialFluidState boundaryFluidState(unsigned globalDofIdx, const int directionId) const
760 {
761 OPM_TIMEBLOCK_LOCAL(boundaryFluidState, Subsystem::Assembly);
762 const auto& bcprop = this->simulator().vanguard().schedule()[this->episodeIndex()].bcprop;
763 if (bcprop.size() > 0) {
764 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
765
766 // index == 0: no boundary conditions for this
767 // global cell and direction
768 if (this->bcindex_(dir)[globalDofIdx] == 0)
769 return initialFluidStates_[globalDofIdx];
770
771 const auto& bc = bcprop[this->bcindex_(dir)[globalDofIdx]];
772 if (bc.bctype == BCType::DIRICHLET )
773 {
774 InitialFluidState fluidState;
775 const int pvtRegionIdx = this->pvtRegionIndex(globalDofIdx);
776 fluidState.setPvtRegionIndex(pvtRegionIdx);
777
778 switch (bc.component) {
779 case BCComponent::OIL:
780 if (!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx))
781 throw std::logic_error("oil is not active and you're trying to add oil BC");
782
783 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 1.0);
784 break;
785 case BCComponent::GAS:
786 if (!FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
787 throw std::logic_error("gas is not active and you're trying to add gas BC");
788
789 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 1.0);
790 break;
791 case BCComponent::WATER:
792 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
793 throw std::logic_error("water is not active and you're trying to add water BC");
794
795 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 1.0);
796 break;
797 case BCComponent::SOLVENT:
798 case BCComponent::POLYMER:
799 case BCComponent::MICR:
800 case BCComponent::OXYG:
801 case BCComponent::UREA:
803 throw std::logic_error("you need to specify a valid component (OIL, WATER or GAS) when DIRICHLET type is set in BC");
804 }
805 fluidState.setTotalSaturation(1.0);
806 double pressure = initialFluidStates_[globalDofIdx].pressure(this->refPressurePhaseIdx_());
807 const auto pressure_input = bc.pressure;
808 if (pressure_input) {
809 pressure = *pressure_input;
810 }
811
812 std::array<Scalar, numPhases> pc = {0};
813 const auto& matParams = this->materialLawParams(globalDofIdx);
814 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
815 Valgrind::CheckDefined(pressure);
816 Valgrind::CheckDefined(pc);
817 for (unsigned activePhaseIdx = 0; activePhaseIdx < FluidSystem::numActivePhases(); ++activePhaseIdx) {
818 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
819 if (Indices::oilEnabled)
820 fluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
821 else if (Indices::gasEnabled)
822 fluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
823 else if (Indices::waterEnabled)
824 //single (water) phase
825 fluidState.setPressure(phaseIdx, pressure);
826 }
827 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
828 double temperature = initialFluidStates_[globalDofIdx].temperature(0); // we only have one temperature
829 const auto temperature_input = bc.temperature;
830 if(temperature_input)
831 temperature = *temperature_input;
832 fluidState.setTemperature(temperature);
833 }
834
835 if constexpr (enableDissolvedGas) {
836 if (FluidSystem::enableDissolvedGas()) {
837 fluidState.setRs(0.0);
838 fluidState.setRv(0.0);
839 }
840 }
841 if constexpr (enableDisgasInWater) {
842 if (FluidSystem::enableDissolvedGasInWater()) {
843 fluidState.setRsw(0.0);
844 }
845 }
846 if constexpr (enableVapwat) {
847 if (FluidSystem::enableVaporizedWater()) {
848 fluidState.setRvw(0.0);
849 }
850 }
851
852 for (unsigned activePhaseIdx = 0; activePhaseIdx < FluidSystem::numActivePhases(); ++activePhaseIdx) {
853 const auto phaseIdx = FluidSystem::activeToCanonicalPhaseIdx(activePhaseIdx);
854
855 const auto& b = FluidSystem::inverseFormationVolumeFactor(fluidState, phaseIdx, pvtRegionIdx);
856 fluidState.setInvB(phaseIdx, b);
857
858 const auto& rho = FluidSystem::density(fluidState, phaseIdx, pvtRegionIdx);
859 fluidState.setDensity(phaseIdx, rho);
860 if constexpr (energyModuleType == EnergyModules::FullyImplicitThermal) {
861 const auto& h = FluidSystem::enthalpy(fluidState, phaseIdx, pvtRegionIdx);
862 fluidState.setEnthalpy(phaseIdx, h);
863 }
864 }
865 fluidState.checkDefined();
866 return fluidState;
867 }
868 }
869 return initialFluidStates_[globalDofIdx];
870 }
871
872
874 { return *eclWriter_; }
875
877 { return *eclWriter_; }
878
883 Scalar maxGasDissolutionFactor(unsigned timeIdx, unsigned globalDofIdx) const
884 {
885 return this->mixControls_.maxGasDissolutionFactor(timeIdx, globalDofIdx,
886 this->episodeIndex(),
887 this->pvtRegionIndex(globalDofIdx));
888 }
889
894 Scalar maxOilVaporizationFactor(unsigned timeIdx, unsigned globalDofIdx) const
895 {
896 return this->mixControls_.maxOilVaporizationFactor(timeIdx, globalDofIdx,
897 this->episodeIndex(),
898 this->pvtRegionIndex(globalDofIdx));
899 }
900
910 {
911 int episodeIdx = this->episodeIndex();
912 return !this->mixControls_.drsdtActive(episodeIdx) &&
913 !this->mixControls_.drvdtActive(episodeIdx) &&
914 this->rockCompPoroMultWc_.empty() &&
915 this->rockCompPoroMult_.empty();
916 }
917
924 template <class Context>
925 void initial(PrimaryVariables& values, const Context& context, unsigned spaceIdx, unsigned timeIdx) const
926 {
927 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
928
929 values.setPvtRegionIndex(pvtRegionIndex(context, spaceIdx, timeIdx));
930 values.assignNaive(initialFluidStates_[globalDofIdx]);
931
932 SolventModule::assignPrimaryVars(values,
933 enableSolvent ? this->solventSaturation_[globalDofIdx] : 0.0,
934 enableSolvent ? this->solventRsw_[globalDofIdx] : 0.0);
935
936 if constexpr (enablePolymer)
937 values[Indices::polymerConcentrationIdx] = this->polymer_.concentration[globalDofIdx];
938
939 if constexpr (enablePolymerMolarWeight)
940 values[Indices::polymerMoleWeightIdx]= this->polymer_.moleWeight[globalDofIdx];
941
942 if constexpr (enableBrine) {
943 if (enableSaltPrecipitation && values.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) {
944 values[Indices::saltConcentrationIdx] = initialFluidStates_[globalDofIdx].saltSaturation();
945 }
946 else {
947 values[Indices::saltConcentrationIdx] = initialFluidStates_[globalDofIdx].saltConcentration();
948 }
949 }
950
951 if constexpr (enableBioeffects) {
952 values[Indices::microbialConcentrationIdx] = this->bioeffects_.microbialConcentration[globalDofIdx];
953 values[Indices::biofilmVolumeFractionIdx]= this->bioeffects_.biofilmVolumeFraction[globalDofIdx];
954 if constexpr (enableMICP) {
955 values[Indices::oxygenConcentrationIdx]= this->bioeffects_.oxygenConcentration[globalDofIdx];
956 values[Indices::ureaConcentrationIdx]= this->bioeffects_.ureaConcentration[globalDofIdx];
957 values[Indices::calciteVolumeFractionIdx]= this->bioeffects_.calciteVolumeFraction[globalDofIdx];
958 }
959 }
960
961 values.checkDefined();
962 }
963
964
965 Scalar drsdtcon(unsigned elemIdx, int episodeIdx) const
966 {
967 return this->mixControls_.drsdtcon(elemIdx, episodeIdx,
968 this->pvtRegionIndex(elemIdx));
969 }
970
971 bool drsdtconIsActive(unsigned elemIdx, int episodeIdx) const
972 {
973 return this->mixControls_.drsdtConvective(episodeIdx, this->pvtRegionIndex(elemIdx));
974 }
975
981 template <class Context>
982 void boundary(BoundaryRateVector& values,
983 const Context& context,
984 unsigned spaceIdx,
985 unsigned timeIdx) const
986 {
987 OPM_TIMEBLOCK_LOCAL(eclProblemBoundary, Subsystem::Assembly);
988 if (!context.intersection(spaceIdx).boundary())
989 return;
990
991 if constexpr (energyModuleType != EnergyModules::FullyImplicitThermal || !enableThermalFluxBoundaries)
992 values.setNoFlow();
993 else {
994 // in the energy case we need to specify a non-trivial boundary condition
995 // because the geothermal gradient needs to be maintained. for this, we
996 // simply assume the initial temperature at the boundary and specify the
997 // thermal flow accordingly. in this context, "thermal flow" means energy
998 // flow due to a temerature gradient while assuming no-flow for mass
999 unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
1000 unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
1001 values.setThermalFlow(context, spaceIdx, timeIdx, this->initialFluidStates_[globalDofIdx] );
1002 }
1003
1004 if (this->nonTrivialBoundaryConditions()) {
1005 unsigned indexInInside = context.intersection(spaceIdx).indexInInside();
1006 unsigned interiorDofIdx = context.interiorScvIndex(spaceIdx, timeIdx);
1007 unsigned globalDofIdx = context.globalSpaceIndex(interiorDofIdx, timeIdx);
1008 unsigned pvtRegionIdx = pvtRegionIndex(context, spaceIdx, timeIdx);
1009 const auto [type, massrate] = this->boundaryCondition(globalDofIdx, indexInInside);
1010 if (type == BCType::THERMAL)
1011 values.setThermalFlow(context, spaceIdx, timeIdx, this->boundaryFluidState(globalDofIdx, indexInInside));
1012 else if (type == BCType::FREE || type == BCType::DIRICHLET)
1013 values.setFreeFlow(context, spaceIdx, timeIdx, this->boundaryFluidState(globalDofIdx, indexInInside));
1014 else if (type == BCType::RATE)
1015 values.setMassRate(massrate, pvtRegionIdx);
1016 }
1017 }
1018
1023 void readSolutionFromOutputModule(const int restart_step, bool fip_init)
1024 {
1025 auto& simulator = this->simulator();
1026 const auto& eclState = simulator.vanguard().eclState();
1027
1028 std::size_t numElems = this->model().numGridDof();
1029 this->initialFluidStates_.resize(numElems);
1030 if constexpr (enableSolvent) {
1031 this->solventSaturation_.resize(numElems, 0.0);
1032 this->solventRsw_.resize(numElems, 0.0);
1033 }
1034
1035 if constexpr (enablePolymer)
1036 this->polymer_.concentration.resize(numElems, 0.0);
1037
1038 if constexpr (enablePolymerMolarWeight) {
1039 const std::string msg {"Support of the RESTART for polymer molecular weight "
1040 "is not implemented yet. The polymer weight value will be "
1041 "zero when RESTART begins"};
1042 OpmLog::warning("NO_POLYMW_RESTART", msg);
1043 this->polymer_.moleWeight.resize(numElems, 0.0);
1044 }
1045
1046 if constexpr (enableBioeffects) {
1047 this->bioeffects_.resize(numElems);
1048 }
1049
1050 // Initialize mixing controls before trying to set any lastRx valuesx
1051 this->mixControls_.init(numElems, restart_step, eclState.runspec().tabdims().getNumPVTTables());
1052
1053 if constexpr (enableBioeffects) {
1054 this->bioeffects_ = this->eclWriter_->outputModule().getBioeffects().getSolution();
1055 }
1056
1057 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1058 auto& elemFluidState = this->initialFluidStates_[elemIdx];
1059 elemFluidState.setPvtRegionIndex(pvtRegionIndex(elemIdx));
1060 this->eclWriter_->outputModule().initHysteresisParams(simulator, elemIdx);
1061 this->eclWriter_->outputModule().assignToFluidState(elemFluidState, elemIdx);
1062
1063 // Note: Function processRestartSaturations_() mutates the
1064 // 'ssol' argument--the value from the restart file--if solvent
1065 // is enabled. Then, store the updated solvent saturation into
1066 // 'solventSaturation_'. Otherwise, just pass a dummy value to
1067 // the function and discard the unchanged result. Do not index
1068 // into 'solventSaturation_' unless solvent is enabled.
1069 {
1070 auto ssol = enableSolvent
1071 ? this->eclWriter_->outputModule().getSolventSaturation(elemIdx)
1072 : Scalar(0);
1073
1074 this->processRestartSaturations_(elemFluidState, ssol);
1075
1076 if constexpr (enableSolvent) {
1077 this->solventSaturation_[elemIdx] = ssol;
1078 this->solventRsw_[elemIdx] = this->eclWriter_->outputModule().getSolventRsw(elemIdx);
1079 }
1080 }
1081
1082 // For CO2STORE and H2STORE we need to set the initial temperature for isothermal simulations
1083 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
1084 bool needTemperature = (eclState.runspec().co2Storage() || eclState.runspec().h2Storage());
1085 if (needTemperature) {
1086 const auto& fp = simulator.vanguard().eclState().fieldProps();
1087 elemFluidState.setTemperature(fp.get_double("TEMPI")[elemIdx]);
1088 }
1089 }
1090
1091 this->mixControls_.updateLastValues(elemIdx, elemFluidState.Rs(), elemFluidState.Rv());
1092
1093 if constexpr (enablePolymer)
1094 this->polymer_.concentration[elemIdx] = this->eclWriter_->outputModule().getPolymerConcentration(elemIdx);
1095 // if we need to restart for polymer molecular weight simulation, we need to add related here
1096 }
1097
1098 const int episodeIdx = this->episodeIndex();
1099 this->mixControls_.updateMaxValues(episodeIdx, simulator.timeStepSize());
1100
1101 // assign the restart solution to the current solution. note that we still need
1102 // to compute real initial solution after this because the initial fluid states
1103 // need to be correct for stuff like boundary conditions.
1104 auto& sol = this->model().solution(/*timeIdx=*/0);
1105 const auto& gridView = this->gridView();
1106 ElementContext elemCtx(simulator);
1107 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
1108 elemCtx.updatePrimaryStencil(elem);
1109 int elemIdx = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
1110 this->initial(sol[elemIdx], elemCtx, /*spaceIdx=*/0, /*timeIdx=*/0);
1111 }
1112
1113 // make sure that the ghost and overlap entities exhibit the correct
1114 // solution. alternatively, this could be done in the loop above by also
1115 // considering non-interior elements. Since the initial() method might not work
1116 // 100% correctly for such elements, let's play safe and explicitly synchronize
1117 // using message passing.
1118 this->model().syncOverlap();
1119
1120 if (fip_init) {
1121 this->updateReferencePorosity_();
1122 this->mixControls_.init(this->model().numGridDof(),
1123 this->episodeIndex(),
1124 eclState.runspec().tabdims().getNumPVTTables());
1125 }
1126 }
1127
1131 Scalar thresholdPressure(unsigned elem1Idx, unsigned elem2Idx) const
1132 { return thresholdPressures_.thresholdPressure(elem1Idx, elem2Idx); }
1133
1135 { return thresholdPressures_; }
1136
1138 { return thresholdPressures_; }
1139
1141 {
1142 return moduleParams_;
1143 }
1144
1145 template<class Serializer>
1146 void serializeOp(Serializer& serializer)
1147 {
1148 serializer(static_cast<FlowProblemType&>(*this));
1149 serializer(mixControls_);
1150 serializer(*eclWriter_);
1151 }
1152
1153protected:
1154 void updateExplicitQuantities_(int episodeIdx, int timeStepSize, const bool first_step_after_restart) override
1155 {
1156 this->updateExplicitQuantities_(first_step_after_restart);
1157
1158 if constexpr (getPropValue<TypeTag, Properties::EnablePolymer>())
1159 updateMaxPolymerAdsorption_();
1160
1161 mixControls_.updateExplicitQuantities(episodeIdx, timeStepSize);
1162 }
1163
1165 {
1166 // we need to update the max polymer adsoption data for all elements
1167 this->updateProperty_("FlowProblemBlackoil::updateMaxPolymerAdsorption_() failed:",
1168 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1169 {
1170 this->updateMaxPolymerAdsorption_(compressedDofIdx,iq);
1171 });
1172 }
1173
1174 bool updateMaxPolymerAdsorption_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1175 {
1176 const Scalar pa = scalarValue(iq.polymerAdsorption());
1177 auto& mpa = this->polymer_.maxAdsorption;
1178 if (mpa[compressedDofIdx] < pa) {
1179 mpa[compressedDofIdx] = pa;
1180 return true;
1181 } else {
1182 return false;
1183 }
1184 }
1185
1187 {
1188 std::vector<Scalar> sumInvB(numPhases, 0.0);
1189 const auto& gridView = this->gridView();
1190 ElementContext elemCtx(this->simulator());
1191 for(const auto& elem: elements(gridView, Dune::Partitions::interior)) {
1192 elemCtx.updatePrimaryStencil(elem);
1193 int elemIdx = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
1194 const auto& dofFluidState = this->initialFluidStates_[elemIdx];
1195 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1196 if (!FluidSystem::phaseIsActive(phaseIdx))
1197 continue;
1198
1199 sumInvB[phaseIdx] += dofFluidState.invB(phaseIdx);
1200 }
1201 }
1202
1203 std::size_t numDof = this->model().numGridDof();
1204 const auto& comm = this->simulator().vanguard().grid().comm();
1205 comm.sum(sumInvB.data(),sumInvB.size());
1206 Scalar numTotalDof = comm.sum(numDof);
1207
1208 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1209 if (!FluidSystem::phaseIsActive(phaseIdx))
1210 continue;
1211
1212 Scalar avgB = numTotalDof / sumInvB[phaseIdx];
1213 const unsigned solventCompIdx = FluidSystem::solventComponentIndex(phaseIdx);
1214 const unsigned activeSolventCompIdx = FluidSystem::canonicalToActiveCompIdx(solventCompIdx);
1215 this->model().setEqWeight(activeSolventCompIdx, avgB);
1216 }
1217 }
1218
1219 // update the parameters needed for DRSDT and DRVDT
1221 {
1222 OPM_TIMEBLOCK(updateCompositionChangeLimits);
1223 // update the "last Rs" values for all elements, including the ones in the ghost
1224 // and overlap regions
1225 int episodeIdx = this->episodeIndex();
1226 std::array<bool,3> active{this->mixControls_.drsdtConvective(episodeIdx),
1227 this->mixControls_.drsdtActive(episodeIdx),
1228 this->mixControls_.drvdtActive(episodeIdx)};
1229 if (!active[0] && !active[1] && !active[2]) {
1230 return false;
1231 }
1232
1233 this->updateProperty_("FlowProblemBlackoil::updateCompositionChangeLimits_()) failed:",
1234 [this,episodeIdx,active](unsigned compressedDofIdx,
1235 const IntensiveQuantities& iq)
1236 {
1237 const DimMatrix& perm = this->intrinsicPermeability(compressedDofIdx);
1238 const Scalar distZ = active[0] ? this->simulator().vanguard().cellThickness(compressedDofIdx) : 0.0;
1239 const int pvtRegionIdx = this->pvtRegionIndex(compressedDofIdx);
1240 this->mixControls_.update(compressedDofIdx,
1241 iq,
1242 episodeIdx,
1243 this->gravity_[dim - 1],
1244 perm[dim - 1][dim - 1],
1245 distZ,
1246 pvtRegionIdx);
1247 }
1248 );
1249
1250 return true;
1251 }
1252
1254 {
1255 // Throw an exception if the grid has LGRs. Refined grid are not supported for restart.
1256 if(this->simulator().vanguard().grid().maxLevel() > 0) {
1257 throw std::invalid_argument("Refined grids are not yet supported for restart ");
1258 }
1259
1260 // Set the start time of the simulation
1261 auto& simulator = this->simulator();
1262 const auto& schedule = simulator.vanguard().schedule();
1263 const auto& eclState = simulator.vanguard().eclState();
1264 const auto& initconfig = eclState.getInitConfig();
1265 const int restart_step = initconfig.getRestartStep();
1266 {
1267 simulator.setTime(schedule.seconds(restart_step));
1268
1269 simulator.startNextEpisode(simulator.startTime() + simulator.time(),
1270 schedule.stepLength(restart_step));
1271 simulator.setEpisodeIndex(restart_step);
1272 }
1273 this->eclWriter_->beginRestart();
1274
1275 Scalar dt = std::min(this->eclWriter_->restartTimeStepSize(), simulator.episodeLength());
1276 simulator.setTimeStepSize(dt);
1277
1278 this->readSolutionFromOutputModule(restart_step, false);
1279
1280 this->eclWriter_->endRestart();
1281 }
1282
1284 {
1285 const auto& simulator = this->simulator();
1286
1287 // initial condition corresponds to hydrostatic conditions.
1288 EquilInitializer<TypeTag> equilInitializer(simulator, *(this->materialLawManager_));
1289
1290 std::size_t numElems = this->model().numGridDof();
1291 this->initialFluidStates_.resize(numElems);
1292 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1293 auto& elemFluidState = this->initialFluidStates_[elemIdx];
1294 elemFluidState.assign(equilInitializer.initialFluidState(elemIdx));
1295 }
1296 }
1297
1299 {
1300 const auto& simulator = this->simulator();
1301 const auto& vanguard = simulator.vanguard();
1302 const auto& eclState = vanguard.eclState();
1303 const auto& fp = eclState.fieldProps();
1304 bool has_swat = fp.has_double("SWAT");
1305 bool has_sgas = fp.has_double("SGAS");
1306 bool has_rs = fp.has_double("RS");
1307 bool has_rsw = fp.has_double("RSW");
1308 bool has_rv = fp.has_double("RV");
1309 bool has_rvw = fp.has_double("RVW");
1310 bool has_pressure = fp.has_double("PRESSURE");
1311 bool has_salt = fp.has_double("SALT");
1312 bool has_saltp = fp.has_double("SALTP");
1313
1314 // make sure all required quantities are enables
1315 if (Indices::numPhases > 1) {
1316 if (FluidSystem::phaseIsActive(waterPhaseIdx) && !has_swat)
1317 throw std::runtime_error("The ECL input file requires the presence of the SWAT keyword if "
1318 "the water phase is active");
1319 if (FluidSystem::phaseIsActive(gasPhaseIdx) && !has_sgas && FluidSystem::phaseIsActive(oilPhaseIdx))
1320 throw std::runtime_error("The ECL input file requires the presence of the SGAS keyword if "
1321 "the gas phase is active");
1322 }
1323 if (!has_pressure)
1324 throw std::runtime_error("The ECL input file requires the presence of the PRESSURE "
1325 "keyword if the model is initialized explicitly");
1326 if (FluidSystem::enableDissolvedGas() && !has_rs)
1327 throw std::runtime_error("The ECL input file requires the RS keyword to be present if"
1328 " dissolved gas is enabled and the model is initialized explicitly");
1329 if (FluidSystem::enableDissolvedGasInWater() && !has_rsw)
1330 OpmLog::warning("The model is initialized explicitly and the RSW keyword is not present in the"
1331 " ECL input file. The RSW values are set equal to 0");
1332 if (FluidSystem::enableVaporizedOil() && !has_rv)
1333 throw std::runtime_error("The ECL input file requires the RV keyword to be present if"
1334 " vaporized oil is enabled and the model is initialized explicitly");
1335 if (FluidSystem::enableVaporizedWater() && !has_rvw)
1336 throw std::runtime_error("The ECL input file requires the RVW keyword to be present if"
1337 " vaporized water is enabled and the model is initialized explicitly");
1338 if (enableBrine && !has_salt)
1339 throw std::runtime_error("The ECL input file requires the SALT keyword to be present if"
1340 " brine is enabled and the model is initialized explicitly");
1341 if (enableSaltPrecipitation && !has_saltp)
1342 throw std::runtime_error("The ECL input file requires the SALTP keyword to be present if"
1343 " salt precipitation is enabled and the model is initialized explicitly");
1344
1345 std::size_t numDof = this->model().numGridDof();
1346
1347 initialFluidStates_.resize(numDof);
1348
1349 std::vector<double> waterSaturationData;
1350 std::vector<double> gasSaturationData;
1351 std::vector<double> pressureData;
1352 std::vector<double> rsData;
1353 std::vector<double> rswData;
1354 std::vector<double> rvData;
1355 std::vector<double> rvwData;
1356 std::vector<double> tempiData;
1357 std::vector<double> saltData;
1358 std::vector<double> saltpData;
1359
1360 if (FluidSystem::phaseIsActive(waterPhaseIdx) && Indices::numPhases > 1)
1361 waterSaturationData = fp.get_double("SWAT");
1362 else
1363 waterSaturationData.resize(numDof);
1364
1365 if (FluidSystem::phaseIsActive(gasPhaseIdx) && FluidSystem::phaseIsActive(oilPhaseIdx))
1366 gasSaturationData = fp.get_double("SGAS");
1367 else
1368 gasSaturationData.resize(numDof);
1369
1370 pressureData = fp.get_double("PRESSURE");
1371 if (FluidSystem::enableDissolvedGas())
1372 rsData = fp.get_double("RS");
1373
1374 if (FluidSystem::enableDissolvedGasInWater() && has_rsw)
1375 rswData = fp.get_double("RSW");
1376
1377 if (FluidSystem::enableVaporizedOil())
1378 rvData = fp.get_double("RV");
1379
1380 if (FluidSystem::enableVaporizedWater())
1381 rvwData = fp.get_double("RVW");
1382
1383 // initial reservoir temperature
1384 tempiData = fp.get_double("TEMPI");
1385
1386 // initial salt concentration data
1387 if constexpr (enableBrine)
1388 saltData = fp.get_double("SALT");
1389
1390 // initial precipitated salt saturation data
1391 if constexpr (enableSaltPrecipitation)
1392 saltpData = fp.get_double("SALTP");
1393
1394 // calculate the initial fluid states
1395 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1396 auto& dofFluidState = initialFluidStates_[dofIdx];
1397
1398 dofFluidState.setPvtRegionIndex(pvtRegionIndex(dofIdx));
1399
1401 // set temperature
1403 if constexpr (energyModuleType != EnergyModules::NoTemperature) {
1404 Scalar temperatureLoc = tempiData[dofIdx];
1405 if (!std::isfinite(temperatureLoc) || temperatureLoc <= 0)
1406 temperatureLoc = FluidSystem::surfaceTemperature;
1407 dofFluidState.setTemperature(temperatureLoc);
1408 }
1409
1411 // set salt concentration
1413 if constexpr (enableBrine)
1414 dofFluidState.setSaltConcentration(saltData[dofIdx]);
1415
1417 // set precipitated salt saturation
1419 if constexpr (enableSaltPrecipitation)
1420 dofFluidState.setSaltSaturation(saltpData[dofIdx]);
1421
1423 // set saturations
1425 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx))
1426 dofFluidState.setSaturation(FluidSystem::waterPhaseIdx,
1427 waterSaturationData[dofIdx]);
1428
1429 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)){
1430 if (!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)){
1431 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
1432 1.0
1433 - waterSaturationData[dofIdx]);
1434 }
1435 else
1436 dofFluidState.setSaturation(FluidSystem::gasPhaseIdx,
1437 gasSaturationData[dofIdx]);
1438 }
1439 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1440 const Scalar soil = 1.0 - waterSaturationData[dofIdx] - gasSaturationData[dofIdx];
1441 if (soil < smallSaturationTolerance_) {
1442 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
1443 }
1444 else {
1445 dofFluidState.setSaturation(FluidSystem::oilPhaseIdx, soil);
1446 }
1447 }
1448
1450 // set phase pressures
1452 Scalar pressure = pressureData[dofIdx]; // oil pressure (or gas pressure for water-gas system or water pressure for single phase)
1453
1454 // this assumes that capillary pressures only depend on the phase saturations
1455 // and possibly on temperature. (this is always the case for ECL problems.)
1456 std::array<Scalar, numPhases> pc = {0};
1457 const auto& matParams = this->materialLawParams(dofIdx);
1458 MaterialLaw::capillaryPressures(pc, matParams, dofFluidState);
1459 Valgrind::CheckDefined(pressure);
1460 Valgrind::CheckDefined(pc);
1461 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1462 if (!FluidSystem::phaseIsActive(phaseIdx))
1463 continue;
1464
1465 if (Indices::oilEnabled)
1466 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[oilPhaseIdx]));
1467 else if (Indices::gasEnabled)
1468 dofFluidState.setPressure(phaseIdx, pressure + (pc[phaseIdx] - pc[gasPhaseIdx]));
1469 else if (Indices::waterEnabled)
1470 //single (water) phase
1471 dofFluidState.setPressure(phaseIdx, pressure);
1472 }
1473
1474 if constexpr (enableDissolvedGas) {
1475 if (FluidSystem::enableDissolvedGas())
1476 dofFluidState.setRs(rsData[dofIdx]);
1477 else if (Indices::gasEnabled && Indices::oilEnabled)
1478 dofFluidState.setRs(0.0);
1479 if (FluidSystem::enableVaporizedOil())
1480 dofFluidState.setRv(rvData[dofIdx]);
1481 else if (Indices::gasEnabled && Indices::oilEnabled)
1482 dofFluidState.setRv(0.0);
1483 }
1484
1485 if constexpr (enableDisgasInWater) {
1486 if (FluidSystem::enableDissolvedGasInWater() && has_rsw)
1487 dofFluidState.setRsw(rswData[dofIdx]);
1488 }
1489
1490 if constexpr (enableVapwat) {
1491 if (FluidSystem::enableVaporizedWater())
1492 dofFluidState.setRvw(rvwData[dofIdx]);
1493 }
1494
1496 // set invB_
1498 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1499 if (!FluidSystem::phaseIsActive(phaseIdx))
1500 continue;
1501
1502 const auto& b = FluidSystem::inverseFormationVolumeFactor(dofFluidState, phaseIdx, pvtRegionIndex(dofIdx));
1503 dofFluidState.setInvB(phaseIdx, b);
1504
1505 const auto& rho = FluidSystem::density(dofFluidState, phaseIdx, pvtRegionIndex(dofIdx));
1506 dofFluidState.setDensity(phaseIdx, rho);
1507
1508 }
1509 }
1510 }
1511
1512
1513 void processRestartSaturations_(InitialFluidState& elemFluidState, Scalar& solventSaturation)
1514 {
1515 // each phase needs to be above certain value to be claimed to be existing
1516 // this is used to recover some RESTART running with the defaulted single-precision format
1517 Scalar sumSaturation = 0.0;
1518 for (std::size_t phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1519 if (FluidSystem::phaseIsActive(phaseIdx)) {
1520 if (elemFluidState.saturation(phaseIdx) < smallSaturationTolerance_)
1521 elemFluidState.setSaturation(phaseIdx, 0.0);
1522
1523 sumSaturation += elemFluidState.saturation(phaseIdx);
1524 }
1525
1526 }
1527 if constexpr (enableSolvent) {
1528 if (solventSaturation < smallSaturationTolerance_)
1529 solventSaturation = 0.0;
1530
1531 sumSaturation += solventSaturation;
1532 }
1533
1534 assert(sumSaturation > 0.0);
1535
1536 for (std::size_t phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
1537 if (FluidSystem::phaseIsActive(phaseIdx)) {
1538 const Scalar saturation = elemFluidState.saturation(phaseIdx) / sumSaturation;
1539 elemFluidState.setSaturation(phaseIdx, saturation);
1540 }
1541 }
1542 if constexpr (enableSolvent) {
1543 solventSaturation = solventSaturation / sumSaturation;
1544 }
1545 }
1546
1548 {
1549 FlowProblemType::readInitialCondition_();
1550
1551 if constexpr (enableSolvent || enablePolymer || enablePolymerMolarWeight || enableBioeffects)
1552 this->readBlackoilExtentionsInitialConditions_(this->model().numGridDof(),
1553 enableSolvent,
1554 enablePolymer,
1555 enablePolymerMolarWeight,
1556 enableBioeffects,
1557 enableMICP);
1558
1559 }
1560
1561 void handleSolventBC(const BCProp::BCFace& bc, RateVector& rate) const override
1562 {
1563 if constexpr (!enableSolvent)
1564 throw std::logic_error("solvent is disabled and you're trying to add solvent to BC");
1565
1566 rate[Indices::solventSaturationIdx] = bc.rate;
1567 }
1568
1569 void handlePolymerBC(const BCProp::BCFace& bc, RateVector& rate) const override
1570 {
1571 if constexpr (!enablePolymer)
1572 throw std::logic_error("polymer is disabled and you're trying to add polymer to BC");
1573
1574 rate[Indices::polymerConcentrationIdx] = bc.rate;
1575 }
1576
1577 void handleMicrBC(const BCProp::BCFace& bc, RateVector& rate) const override
1578 {
1579 if constexpr (!enableMICP)
1580 throw std::logic_error("MICP is disabled and you're trying to add microbes to BC");
1581
1582 rate[Indices::microbialConcentrationIdx] = bc.rate;
1583 }
1584
1585 void handleOxygBC(const BCProp::BCFace& bc, RateVector& rate) const override
1586 {
1587 if constexpr (!enableMICP)
1588 throw std::logic_error("MICP is disabled and you're trying to add oxygen to BC");
1589
1590 rate[Indices::oxygenConcentrationIdx] = bc.rate;
1591 }
1592
1593 void handleUreaBC(const BCProp::BCFace& bc, RateVector& rate) const override
1594 {
1595 if constexpr (!enableMICP)
1596 throw std::logic_error("MICP is disabled and you're trying to add urea to BC");
1597
1598 rate[Indices::ureaConcentrationIdx] = bc.rate;
1599 // since the urea concentration can be much larger than 1, then we apply a scaling factor
1600 rate[Indices::ureaConcentrationIdx] *= getPropValue<TypeTag, Properties::BlackOilUreaScalingFactor>();
1601 }
1602
1603 void updateExplicitQuantities_(const bool first_step_after_restart)
1604 {
1605 OPM_TIMEBLOCK(updateExplicitQuantities);
1606 const bool invalidateFromMaxWaterSat = this->updateMaxWaterSaturation_();
1607 const bool invalidateFromMinPressure = this->updateMinPressure_();
1608
1609 // update hysteresis and max oil saturation used in vappars
1610 const bool invalidateFromHyst = this->updateHysteresis_();
1611 const bool invalidateFromMaxOilSat = this->updateMaxOilSaturation_();
1612
1613 // deal with DRSDT and DRVDT
1614 const bool invalidateDRDT = !first_step_after_restart && this->updateCompositionChangeLimits_();
1615
1616 // the derivatives may have changed
1617 const bool invalidateIntensiveQuantities
1618 = invalidateFromMaxWaterSat || invalidateFromMinPressure || invalidateFromHyst || invalidateFromMaxOilSat || invalidateDRDT;
1619 if (invalidateIntensiveQuantities) {
1620 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
1621 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
1622 }
1623
1624 this->updateRockCompTransMultVal_();
1625 }
1626
1628 {
1629 if (const auto nph = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
1630 + FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)
1631 + FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
1632 nph < 2)
1633 {
1634 // Single phase runs don't need saturation functions and there's
1635 // nothing to do here. Return 'true' to tell caller that the
1636 // consistency requirements are Met.
1637 return true;
1638 }
1639
1640 const auto numSamplePoints = static_cast<std::size_t>
1641 (Parameters::Get<Parameters::NumSatfuncConsistencySamplePoints>());
1642
1643 auto sfuncConsistencyChecks =
1645 numSamplePoints, this->simulator().vanguard().eclState(),
1646 [&cmap = this->simulator().vanguard().cartesianIndexMapper()](const int elemIdx)
1647 { return cmap.cartesianIndex(elemIdx); }
1648 };
1649
1650 const auto ioRank = 0;
1651 const auto isIoRank = this->simulator().vanguard()
1652 .grid().comm().rank() == ioRank;
1653
1654 // Note: Run saturation function consistency checks on main grid
1655 // only (i.e., levelGridView(0)). These checks are not supported
1656 // for LGRs at this time.
1657 sfuncConsistencyChecks.collectFailuresTo(ioRank)
1658 .run(this->simulator().vanguard().grid().levelGridView(0),
1659 [&vg = this->simulator().vanguard(),
1660 &emap = this->simulator().model().elementMapper()]
1661 (const auto& elem)
1662 { return vg.gridIdxToEquilGridIdx(emap.index(elem)); });
1663
1664 using ViolationLevel = typename Satfunc::PhaseChecks::
1665 SatfuncConsistencyCheckManager<Scalar>::ViolationLevel;
1666
1667 auto reportFailures = [&sfuncConsistencyChecks]
1668 (const ViolationLevel level)
1669 {
1670 sfuncConsistencyChecks.reportFailures
1671 (level, [](std::string_view record)
1672 { OpmLog::info(std::string { record }); });
1673 };
1674
1675 if (sfuncConsistencyChecks.anyFailedStandardChecks()) {
1676 if (isIoRank) {
1677 OpmLog::warning("Saturation Function "
1678 "End-point Consistency Problems");
1679
1680 reportFailures(ViolationLevel::Standard);
1681 }
1682 }
1683
1684 if (sfuncConsistencyChecks.anyFailedCriticalChecks()) {
1685 if (isIoRank) {
1686 OpmLog::error("Saturation Function "
1687 "End-point Consistency Failures");
1688
1689 reportFailures(ViolationLevel::Critical);
1690 }
1691
1692 // There are "critical" check failures. Report that consistency
1693 // requirements are not Met.
1694 return false;
1695 }
1696
1697 // If we get here then there are no critical failures. Report
1698 // Met = true, i.e., that the consistency requirements ARE met.
1699 return true;
1700 }
1701
1703
1704 std::vector<InitialFluidState> initialFluidStates_;
1705
1707 std::unique_ptr<EclWriterType> eclWriter_;
1708
1709 const Scalar smallSaturationTolerance_ = 1.e-6;
1710#if HAVE_DAMARIS
1711 bool enableDamarisOutput_ = false ;
1712 std::unique_ptr<DamarisWriterType> damarisWriter_;
1713#endif
1715
1717
1719
1721
1722private:
1733 bool episodeWillBeOver() const override
1734 {
1735 const auto currTime = this->simulator().time()
1736 + this->simulator().timeStepSize();
1737
1738 const auto nextReportStep =
1739 this->simulator().vanguard().schedule()
1740 .seconds(this->simulator().episodeIndex() + 1);
1741
1742 const auto isSubStep = (nextReportStep - currTime)
1743 > (2 * std::numeric_limits<float>::epsilon()) * nextReportStep;
1744
1745 return !isSubStep;
1746 }
1747};
1748
1749} // namespace Opm
1750
1751#endif // OPM_FLOW_PROBLEM_BLACK_HPP
Classes required for dynamic convective mixing.
Class handling Action support in simulator.
Definition: ActionHandler.hpp:52
static void setParams(BlackOilBioeffectsParams< Scalar > &&params)
Set parameters.
Definition: blackoilbioeffectsmodules.hh:133
static void setParams(BlackOilBrineParams< Scalar > &&params)
Set parameters.
Definition: blackoilbrinemodules.hh:90
static void setParams(BlackOilExtboParams< Scalar > &&params)
Set parameters.
Definition: blackoilextbomodules.hh:89
static void setParams(BlackOilFoamParams< Scalar > &&params)
Set parameters.
Definition: blackoilfoammodules.hh:89
Hybrid Newton solver extension for the black-oil model.
Definition: HybridNewton.hpp:59
void tryApplyHybridNewton()
Attempt to apply the Hybrid Newton correction at the current timestep.
Definition: HybridNewton.hpp:99
static void setParams(BlackOilPolymerParams< Scalar > &&params)
Set parameters.
Definition: blackoilpolymermodules.hh:96
static void setParams(BlackOilSolventParams< Scalar > &&params)
Set parameters.
Definition: blackoilsolventmodules.hh:101
Collects necessary output values and pass it to opm-common's ECL output.
Definition: EclWriter.hpp:115
static void registerParameters()
Definition: EclWriter.hpp:138
Computes the initial condition based on the EQUIL keyword from ECL.
Definition: EquilInitializer.hpp:59
BlackOilFluidState< Scalar, FluidSystem, energyModuleType !=EnergyModules::NoTemperature, energyModuleType==EnergyModules::FullyImplicitThermal, enableDissolution, enableVapwat, enableBrine, enableSaltPrecipitation, enableDisgasInWater, Indices::numPhases > ScalarFluidState
Definition: EquilInitializer.hpp:98
const ScalarFluidState & initialFluidState(unsigned elemIdx) const
Return the initial thermodynamic state which should be used as the initial condition.
Definition: EquilInitializer.hpp:198
void readRockParameters_(const std::vector< Scalar > &cellCenterDepths, std::function< std::array< int, 3 >(const unsigned)> ijkIndex)
Definition: FlowGenericProblem_impl.hpp:137
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblemBlackoil.hpp:85
HybridNewton hybridNewton_
Definition: FlowProblemBlackoil.hpp:1720
void updateExplicitQuantities_(int episodeIdx, int timeStepSize, const bool first_step_after_restart) override
Definition: FlowProblemBlackoil.hpp:1154
bool updateMaxPolymerAdsorption_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblemBlackoil.hpp:1174
void handlePolymerBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1569
void writeOutput(const bool verbose) override
Write the requested quantities of the current solution into the output files.
Definition: FlowProblemBlackoil.hpp:559
void readInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1547
void handleOxygBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1585
void readEquilInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1283
void handleSolventBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1561
Scalar maxGasDissolutionFactor(unsigned timeIdx, unsigned globalDofIdx) const
Returns the maximum value of the gas dissolution factor at the current time for a given degree of fre...
Definition: FlowProblemBlackoil.hpp:883
const std::vector< InitialFluidState > & initialFluidStates() const
Definition: FlowProblemBlackoil.hpp:747
void processRestartSaturations_(InitialFluidState &elemFluidState, Scalar &solventSaturation)
Definition: FlowProblemBlackoil.hpp:1513
std::vector< InitialFluidState > & initialFluidStates()
Definition: FlowProblemBlackoil.hpp:744
FlowProblemBlackoil(Simulator &simulator)
Definition: FlowProblemBlackoil.hpp:184
bool enableEclOutput_
Definition: FlowProblemBlackoil.hpp:1706
Scalar drsdtcon(unsigned elemIdx, int episodeIdx) const
Definition: FlowProblemBlackoil.hpp:965
void endStepApplyAction()
Definition: FlowProblemBlackoil.hpp:480
bool drsdtconIsActive(unsigned elemIdx, int episodeIdx) const
Definition: FlowProblemBlackoil.hpp:971
Scalar maxOilVaporizationFactor(unsigned timeIdx, unsigned globalDofIdx) const
Returns the maximum value of the oil vaporization factor at the current time for a given degree of fr...
Definition: FlowProblemBlackoil.hpp:894
std::vector< InitialFluidState > initialFluidStates_
Definition: FlowProblemBlackoil.hpp:1704
void endTimeStep() override
Called by the simulator after each time integration.
Definition: FlowProblemBlackoil.hpp:474
void updateMaxPolymerAdsorption_()
Definition: FlowProblemBlackoil.hpp:1164
const InitialFluidState & initialFluidState(unsigned globalDofIdx) const
Definition: FlowProblemBlackoil.hpp:741
void endEpisode() override
Called by the simulator after the end of an episode.
Definition: FlowProblemBlackoil.hpp:527
void setSubStepReport(const SimulatorReportSingle &report)
Definition: FlowProblemBlackoil.hpp:753
void initial(PrimaryVariables &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the initial value for a control volume.
Definition: FlowProblemBlackoil.hpp:925
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition: FlowProblemBlackoil.hpp:283
void handleMicrBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1577
const FlowThresholdPressure< TypeTag > & thresholdPressure() const
Definition: FlowProblemBlackoil.hpp:1134
void finalizeOutput()
Definition: FlowProblemBlackoil.hpp:580
void boundary(BoundaryRateVector &values, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the boundary conditions for a boundary segment.
Definition: FlowProblemBlackoil.hpp:982
InitialFluidState boundaryFluidState(unsigned globalDofIdx, const int directionId) const
Definition: FlowProblemBlackoil.hpp:759
std::unique_ptr< EclWriterType > eclWriter_
Definition: FlowProblemBlackoil.hpp:1707
void initialSolutionApplied() override
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblemBlackoil.hpp:592
const ModuleParams & moduleParams() const
Definition: FlowProblemBlackoil.hpp:1140
void readEclRestartSolution_()
Definition: FlowProblemBlackoil.hpp:1253
FlowThresholdPressure< TypeTag > & thresholdPressure()
Definition: FlowProblemBlackoil.hpp:1137
FlowThresholdPressure< TypeTag > thresholdPressures_
Definition: FlowProblemBlackoil.hpp:1702
void readExplicitInitialCondition_() override
Definition: FlowProblemBlackoil.hpp:1298
void beginEpisode() override
Called by the simulator before an episode begins.
Definition: FlowProblemBlackoil.hpp:243
bool recycleFirstIterationStorage() const
Return if the storage term of the first iteration is identical to the storage term for the solution o...
Definition: FlowProblemBlackoil.hpp:909
LhsEval permFactTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Calculate the transmissibility multiplier due to porosity reduction.
Definition: FlowProblemBlackoil.hpp:721
void serializeOp(Serializer &serializer)
Definition: FlowProblemBlackoil.hpp:1146
MixingRateControls< FluidSystem > mixControls_
Definition: FlowProblemBlackoil.hpp:1714
void writeReports(const SimulatorTimer &timer)
Definition: FlowProblemBlackoil.hpp:547
ModuleParams moduleParams_
Definition: FlowProblemBlackoil.hpp:1718
const EclWriterType & eclWriter() const
Definition: FlowProblemBlackoil.hpp:873
void setSimulationReport(const SimulatorReport &report)
Definition: FlowProblemBlackoil.hpp:756
void addToSourceDense(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const override
Definition: FlowProblemBlackoil.hpp:614
void computeAndSetEqWeights_()
Definition: FlowProblemBlackoil.hpp:1186
void beginTimeStep() override
Called by the simulator before each time integration.
Definition: FlowProblemBlackoil.hpp:274
void updateExplicitQuantities_(const bool first_step_after_restart)
Definition: FlowProblemBlackoil.hpp:1603
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblemBlackoil.hpp:170
ActionHandler< Scalar, IndexTraits > actionHandler_
Definition: FlowProblemBlackoil.hpp:1716
void readSolutionFromOutputModule(const int restart_step, bool fip_init)
Read simulator solution state from the outputmodule (used with restart)
Definition: FlowProblemBlackoil.hpp:1023
const EclipseIO & eclIO() const
Definition: FlowProblemBlackoil.hpp:750
Scalar thresholdPressure(unsigned elem1Idx, unsigned elem2Idx) const
Definition: FlowProblemBlackoil.hpp:1131
void handleUreaBC(const BCProp::BCFace &bc, RateVector &rate) const override
Definition: FlowProblemBlackoil.hpp:1593
EclWriterType & eclWriter()
Definition: FlowProblemBlackoil.hpp:876
bool satfuncConsistencyRequirementsMet() const
Definition: FlowProblemBlackoil.hpp:1627
bool updateCompositionChangeLimits_()
Definition: FlowProblemBlackoil.hpp:1220
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblem.hpp:95
@ enableSolvent
Definition: FlowProblem.hpp:134
@ enableBioeffects
Definition: FlowProblem.hpp:120
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:495
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:845
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:678
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:108
@ oilPhaseIdx
Definition: FlowProblem.hpp:138
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: FlowProblem.hpp:102
GetPropType< TypeTag, Properties::EqVector > EqVector
Definition: FlowProblem.hpp:107
@ enablePolymerMolarWeight
Definition: FlowProblem.hpp:132
@ gasCompIdx
Definition: FlowProblem.hpp:143
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:151
@ numComponents
Definition: FlowProblem.hpp:118
GlobalEqVector drift_
Definition: FlowProblem.hpp:1776
@ enableExtbo
Definition: FlowProblem.hpp:128
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:148
@ dim
Definition: FlowProblem.hpp:112
@ enableMICP
Definition: FlowProblem.hpp:130
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: FlowProblem.hpp:165
@ numPhases
Definition: FlowProblem.hpp:117
int episodeIndex() const
Definition: FlowProblem.hpp:294
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:109
@ numEq
Definition: FlowProblem.hpp:116
@ oilCompIdx
Definition: FlowProblem.hpp:144
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: FlowProblem.hpp:106
@ enablePolymer
Definition: FlowProblem.hpp:131
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: FlowProblem.hpp:149
@ enableConvectiveMixing
Definition: FlowProblem.hpp:122
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:135
TracerModel tracerModel_
Definition: FlowProblem.hpp:1782
WellModel wellModel_
Definition: FlowProblem.hpp:1778
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:302
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:133
virtual void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:361
@ enableExperiments
Definition: FlowProblem.hpp:127
@ enableFoam
Definition: FlowProblem.hpp:129
void readThermalParameters_()
Definition: FlowProblem.hpp:1421
@ gasPhaseIdx
Definition: FlowProblem.hpp:137
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:160
@ enableDiffusion
Definition: FlowProblem.hpp:123
@ enableBrine
Definition: FlowProblem.hpp:121
TemperatureModel temperatureModel_
Definition: FlowProblem.hpp:1783
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowProblem.hpp:103
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblem.hpp:182
void updatePffDofData_()
Definition: FlowProblem.hpp:1569
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:147
@ waterCompIdx
Definition: FlowProblem.hpp:145
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1600
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1771
static constexpr EnergyModules energyModuleType
Definition: FlowProblem.hpp:125
@ waterPhaseIdx
Definition: FlowProblem.hpp:139
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: FlowProblem.hpp:105
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: FlowProblem.hpp:157
@ dimWorld
Definition: FlowProblem.hpp:113
@ enableDispersion
Definition: FlowProblem.hpp:124
void readMaterialParameters_()
Definition: FlowProblem.hpp:1381
This class calculates the threshold pressure for grid faces according to the Eclipse Reference Manual...
Definition: FlowThresholdPressure.hpp:59
Class handling mixing rate controls for a FlowProblemBlackoil.
Definition: MixingRateControls.hpp:46
Definition: SatfuncConsistencyCheckManager.hpp:58
SatfuncConsistencyCheckManager & collectFailuresTo(const int root)
Definition: SatfuncConsistencyCheckManager.hpp:99
void run(const GridView &gv, GetCellIndex &&getCellIndex)
Definition: SatfuncConsistencyCheckManager.hpp:128
Definition: SimulatorTimer.hpp:39
VTK output module for the tracer model's parameters.
Definition: VtkTracerModule.hpp:58
static void registerParameters()
Register all run-time parameters for the tracer VTK output module.
Definition: VtkTracerModule.hpp:84
@ NONE
Definition: DeferredLogger.hpp:46
static constexpr int dim
Definition: structuredgridvanguard.hh:68
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
Struct holding the parameters for the BlackOilBioeffectsModule class.
Definition: blackoilbioeffectsparams.hpp:42
Struct holding the parameters for the BlackoilBrineModule class.
Definition: blackoilbrineparams.hpp:42
Struct holding the parameters for the BlackoilExtboModule class.
Definition: blackoilextboparams.hpp:47
Struct holding the parameters for the BlackoilFoamModule class.
Definition: blackoilfoamparams.hpp:44
Struct holding the parameters for the BlackOilPolymerModule class.
Definition: blackoilpolymerparams.hpp:43
Struct holding the parameters for the BlackOilSolventModule class.
Definition: blackoilsolventparams.hpp:47
Definition: SimulatorReport.hpp:122
A struct for returning timing data from a simulator to its caller.
Definition: SimulatorReport.hpp:34