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