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