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