FlowProblem.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
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 2 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20
21 Consult the COPYING file in the top-level source directory of this
22 module for the precise wording of the license and the list of
23 copyright holders.
24*/
30#ifndef OPM_FLOW_PROBLEM_HPP
31#define OPM_FLOW_PROBLEM_HPP
32
33#include <dune/common/version.hh>
34#include <dune/common/fvector.hh>
35#include <dune/common/fmatrix.hh>
36
37#include <opm/common/utility/TimeService.hpp>
38
39#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
40#include <opm/input/eclipse/Parser/ParserKeywords/E.hpp>
41#include <opm/input/eclipse/Schedule/Schedule.hpp>
42
43#include <opm/material/common/ConditionalStorage.hpp>
44#include <opm/material/common/Valgrind.hpp>
45#include <opm/material/densead/Evaluation.hpp>
46#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
47#include <opm/material/thermal/EclThermalLawManager.hpp>
48
52
53#include <opm/output/eclipse/EclipseIO.hpp>
54
60// TODO: maybe we can name it FlowProblemProperties.hpp
66
70
71#include <opm/utility/CopyablePtr.hpp>
72
73#include <algorithm>
74#include <functional>
75#include <set>
76#include <string>
77#include <vector>
78
79namespace Opm {
80
87template <class TypeTag>
88class FlowProblem : public GetPropType<TypeTag, Properties::BaseProblem>
89 , public FlowGenericProblem<GetPropType<TypeTag, Properties::GridView>,
90 GetPropType<TypeTag, Properties::FluidSystem>>
91{
92protected:
97
105
106 // Grid and world dimension
107 enum { dim = GridView::dimension };
108 enum { dimWorld = GridView::dimensionworld };
109
110 // copy some indices for convenience
111 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
112 enum { numPhases = FluidSystem::numPhases };
113 enum { numComponents = FluidSystem::numComponents };
114
115 enum { enableConvectiveMixing = getPropValue<TypeTag, Properties::EnableConvectiveMixing>() };
116 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
117 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
118 enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
119 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
120 enum { enableExperiments = getPropValue<TypeTag, Properties::EnableExperiments>() };
121 enum { enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>() };
122 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
123 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
124 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
125 enum { enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>() };
126 enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
127 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
128 enum { enableTemperature = getPropValue<TypeTag, Properties::EnableTemperature>() };
129 enum { enableThermalFluxBoundaries = getPropValue<TypeTag, Properties::EnableThermalFluxBoundaries>() };
130
131 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
132 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
133 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
134
135 // TODO: later, gasCompIdx, oilCompIdx and waterCompIdx should go to the FlowProblemBlackoil in the future
136 // we do not want them in the compositional setting
137 enum { gasCompIdx = FluidSystem::gasCompIdx };
138 enum { oilCompIdx = FluidSystem::oilCompIdx };
139 enum { waterCompIdx = FluidSystem::waterCompIdx };
140
144 using Element = typename GridView::template Codim<0>::Entity;
148 using MaterialLawParams = typename EclMaterialLawManager::MaterialLawParams;
149 using SolidEnergyLawParams = typename EclThermalLawManager::SolidEnergyLawParams;
150 using ThermalConductionLawParams = typename EclThermalLawManager::ThermalConductionLawParams;
158
159 using Toolbox = MathToolbox<Evaluation>;
160 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
161
162
164 using DirectionalMobilityPtr = Utility::CopyablePtr<DirectionalMobility<TypeTag, Evaluation>>;
165
166public:
173 using BaseType::porosity;
174
178 static void registerParameters()
179 {
180 ParentType::registerParameters();
181
182 registerFlowProblemParameters<Scalar>();
183 }
184
185
189 static int handlePositionalParameter(std::function<void(const std::string&,
190 const std::string&)> addKey,
191 std::set<std::string>& seenParams,
192 std::string& errorMsg,
193 int,
194 const char** argv,
195 int paramIdx,
196 int)
197 {
198 return eclPositionalParameter(addKey,
199 seenParams,
200 errorMsg,
201 argv,
202 paramIdx);
203 }
204
208 explicit FlowProblem(Simulator& simulator)
209 : ParentType(simulator)
210 , BaseType(simulator.vanguard().eclState(),
211 simulator.vanguard().schedule(),
212 simulator.vanguard().gridView())
213 , transmissibilities_(simulator.vanguard().eclState(),
214 simulator.vanguard().gridView(),
215 simulator.vanguard().cartesianIndexMapper(),
216 simulator.vanguard().grid(),
217 simulator.vanguard().cellCentroids(),
221 , wellModel_(simulator)
222 , aquiferModel_(simulator)
223 , pffDofData_(simulator.gridView(), this->elementMapper())
224 , tracerModel_(simulator)
225 {
226 const auto& vanguard = simulator.vanguard();
227
228 enableDriftCompensation_ = Parameters::Get<Parameters::EnableDriftCompensation>();
229
230 enableVtkOutput_ = Parameters::Get<Parameters::EnableVtkOutput>();
231
232 this->enableTuning_ = Parameters::Get<Parameters::EnableTuning>();
233 this->initialTimeStepSize_ = Parameters::Get<Parameters::InitialTimeStepSize<Scalar>>();
234 this->maxTimeStepAfterWellEvent_ = Parameters::Get<Parameters::TimeStepAfterEventInDays<Scalar>>() * 24 * 60 * 60;
235
236 // The value N for this parameter is defined in the following order of presedence:
237 // 1. Command line value (--num-pressure-points-equil=N)
238 // 2. EQLDIMS item 2
239 // Default value is defined in opm-common/src/opm/input/eclipse/share/keywords/000_Eclipse100/E/EQLDIMS
240 if (Parameters::IsSet<Parameters::NumPressurePointsEquil>())
241 {
242 this->numPressurePointsEquil_ = Parameters::Get<Parameters::NumPressurePointsEquil>();
243 } else {
244 this->numPressurePointsEquil_ = simulator.vanguard().eclState().getTableManager().getEqldims().getNumDepthNodesP();
245 }
246
247 explicitRockCompaction_ = Parameters::Get<Parameters::ExplicitRockCompaction>();
248
249
250 RelpermDiagnostics relpermDiagnostics;
251 relpermDiagnostics.diagnosis(vanguard.eclState(), vanguard.cartesianIndexMapper());
252 }
253
254 virtual ~FlowProblem() = default;
255
256 void prefetch(const Element& elem) const
257 { pffDofData_.prefetch(elem); }
258
270 template <class Restarter>
271 void deserialize(Restarter& res)
272 {
273 // reload the current episode/report step from the deck
274 this->beginEpisode();
275
276 // deserialize the wells
277 wellModel_.deserialize(res);
278
279 // deserialize the aquifer
280 aquiferModel_.deserialize(res);
281 }
282
289 template <class Restarter>
290 void serialize(Restarter& res)
291 {
292 wellModel_.serialize(res);
293
294 aquiferModel_.serialize(res);
295 }
296
297 int episodeIndex() const
298 {
299 return std::max(this->simulator().episodeIndex(), 0);
300 }
301
305 virtual void beginEpisode()
306 {
307 OPM_TIMEBLOCK(beginEpisode);
308 // Proceed to the next report step
309 auto& simulator = this->simulator();
310 int episodeIdx = simulator.episodeIndex();
311 auto& eclState = simulator.vanguard().eclState();
312 const auto& schedule = simulator.vanguard().schedule();
313 const auto& events = schedule[episodeIdx].events();
314
315 if (episodeIdx >= 0 && events.hasEvent(ScheduleEvents::GEO_MODIFIER)) {
316 // bring the contents of the keywords to the current state of the SCHEDULE
317 // section.
318 //
319 // TODO (?): make grid topology changes possible (depending on what exactly
320 // has changed, the grid may need be re-created which has some serious
321 // implications on e.g., the solution of the simulation.)
322 const auto& miniDeck = schedule[episodeIdx].geo_keywords();
323 const auto& cc = simulator.vanguard().grid().comm();
324 eclState.apply_schedule_keywords( miniDeck );
325 eclBroadcast(cc, eclState.getTransMult() );
326
327 // Re-ordering in case of ALUGrid
328 std::function<unsigned int(unsigned int)> equilGridToGrid = [&simulator](unsigned int i) {
329 return simulator.vanguard().gridEquilIdxToGridIdx(i);
330 };
331
332 // re-compute all quantities which may possibly be affected.
333 using TransUpdateQuantities = typename Vanguard::TransmissibilityType::TransUpdateQuantities;
334 transmissibilities_.update(true, TransUpdateQuantities::All, equilGridToGrid);
335 this->referencePorosity_[1] = this->referencePorosity_[0];
338 this->model().linearizer().updateDiscretizationParameters();
339 }
340
341 bool tuningEvent = this->beginEpisode_(enableExperiments, this->episodeIndex());
342
343 // set up the wells for the next episode.
344 wellModel_.beginEpisode();
345
346 // set up the aquifers for the next episode.
347 aquiferModel_.beginEpisode();
348
349 // set the size of the initial time step of the episode
350 Scalar dt = limitNextTimeStepSize_(simulator.episodeLength());
351 // negative value of initialTimeStepSize_ indicates no active limit from TSINIT or NEXTSTEP
352 if ( (episodeIdx == 0 || tuningEvent) && this->initialTimeStepSize_ > 0)
353 // allow the size of the initial time step to be set via an external parameter
354 // if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
355 dt = std::min(dt, this->initialTimeStepSize_);
356 simulator.setTimeStepSize(dt);
357 }
358
363 {
364 OPM_TIMEBLOCK(beginTimeStep);
365 const int episodeIdx = this->episodeIndex();
366 const int timeStepSize = this->simulator().timeStepSize();
367
369 episodeIdx,
370 this->simulator().timeStepIndex(),
371 this->simulator().startTime(),
372 this->simulator().time(),
373 timeStepSize,
374 this->simulator().endTime());
375
376 // update maximum water saturation and minimum pressure
377 // used when ROCKCOMP is activated
378 // Do not update max RS first step after a restart
379 this->updateExplicitQuantities_(episodeIdx, timeStepSize, first_step_ && (episodeIdx > 0));
380 first_step_ = false;
381
383 this->model().linearizer().updateBoundaryConditionData();
384 }
385
386 wellModel_.beginTimeStep();
387 aquiferModel_.beginTimeStep();
388 tracerModel_.beginTimeStep();
389
390 }
391
396 {
397 OPM_TIMEBLOCK(beginIteration);
398 wellModel_.beginIteration();
399 aquiferModel_.beginIteration();
400 }
401
406 {
407 OPM_TIMEBLOCK(endIteration);
408 wellModel_.endIteration();
409 aquiferModel_.endIteration();
410 }
411
415 virtual void endTimeStep()
416 {
417 OPM_TIMEBLOCK(endTimeStep);
418
419#ifndef NDEBUG
420 if constexpr (getPropValue<TypeTag, Properties::EnableDebuggingChecks>()) {
421 // in debug mode, we don't care about performance, so we check
422 // if the model does the right thing (i.e., the mass change
423 // inside the whole reservoir must be equivalent to the fluxes
424 // over the grid's boundaries plus the source rates specified by
425 // the problem).
426 const int rank = this->simulator().gridView().comm().rank();
427 if (rank == 0) {
428 std::cout << "checking conservativeness of solution\n";
429 }
430
431 this->model().checkConservativeness(/*tolerance=*/-1, /*verbose=*/true);
432 if (rank == 0) {
433 std::cout << "solution is sufficiently conservative\n";
434 }
435 }
436#endif // NDEBUG
437
438 auto& simulator = this->simulator();
439 simulator.setTimeStepIndex(simulator.timeStepIndex()+1);
440
441 this->wellModel_.endTimeStep();
442 this->aquiferModel_.endTimeStep();
443 this->tracerModel_.endTimeStep();
444
445 // Compute flux for output
446 this->model().linearizer().updateFlowsInfo();
447
448 if (this->enableDriftCompensation_) {
449 OPM_TIMEBLOCK(driftCompansation);
450
451 const auto& residual = this->model().linearizer().residual();
452
453 for (unsigned globalDofIdx = 0; globalDofIdx < residual.size(); globalDofIdx ++) {
454 this->drift_[globalDofIdx] = residual[globalDofIdx] * simulator.timeStepSize();
455
456 if constexpr (getPropValue<TypeTag, Properties::UseVolumetricResidual>()) {
457 this->drift_[globalDofIdx] *= this->model().dofTotalVolume(globalDofIdx);
458 }
459 }
460 }
461 }
462
466 virtual void endEpisode()
467 {
468 const int episodeIdx = this->episodeIndex();
469
470 this->wellModel_.endEpisode();
471 this->aquiferModel_.endEpisode();
472
473 const auto& schedule = this->simulator().vanguard().schedule();
474
475 // End simulation when completed.
476 if (episodeIdx + 1 >= static_cast<int>(schedule.size()) - 1) {
477 this->simulator().setFinished(true);
478 return;
479 }
480
481 // Otherwise, start next episode (report step).
482 this->simulator().startNextEpisode(schedule.stepLength(episodeIdx + 1));
483 }
484
489 void writeOutput(bool verbose = true)
490 {
491 OPM_TIMEBLOCK(problemWriteOutput);
492 // use the generic code to prepare the output fields and to
493 // write the desired VTK files.
494 if (Parameters::Get<Parameters::EnableWriteAllSolutions>() ||
495 this->simulator().episodeWillBeOver()) {
496 ParentType::writeOutput(verbose);
497 }
498 }
499
503 template <class Context>
504 const DimMatrix& intrinsicPermeability(const Context& context,
505 unsigned spaceIdx,
506 unsigned timeIdx) const
507 {
508 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
509 return transmissibilities_.permeability(globalSpaceIdx);
510 }
511
518 const DimMatrix& intrinsicPermeability(unsigned globalElemIdx) const
519 { return transmissibilities_.permeability(globalElemIdx); }
520
524 template <class Context>
525 Scalar transmissibility(const Context& context,
526 [[maybe_unused]] unsigned fromDofLocalIdx,
527 unsigned toDofLocalIdx) const
528 {
529 assert(fromDofLocalIdx == 0);
530 return pffDofData_.get(context.element(), toDofLocalIdx).transmissibility;
531 }
532
536 Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
537 {
538 return transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
539 }
540
544 template <class Context>
545 Scalar diffusivity(const Context& context,
546 [[maybe_unused]] unsigned fromDofLocalIdx,
547 unsigned toDofLocalIdx) const
548 {
549 assert(fromDofLocalIdx == 0);
550 return *pffDofData_.get(context.element(), toDofLocalIdx).diffusivity;
551 }
552
556 Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
557 return transmissibilities_.diffusivity(globalCellIn, globalCellOut);
558 }
559
563 Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
564 return transmissibilities_.dispersivity(globalCellIn, globalCellOut);
565 }
566
570 Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx,
571 const unsigned boundaryFaceIdx) const
572 {
573 return transmissibilities_.thermalTransmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
574 }
575
576
577
578
582 template <class Context>
583 Scalar transmissibilityBoundary(const Context& elemCtx,
584 unsigned boundaryFaceIdx) const
585 {
586 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
587 return transmissibilities_.transmissibilityBoundary(elemIdx, boundaryFaceIdx);
588 }
589
593 Scalar transmissibilityBoundary(const unsigned globalSpaceIdx,
594 const unsigned boundaryFaceIdx) const
595 {
596 return transmissibilities_.transmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
597 }
598
599
603 Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn,
604 const unsigned globalSpaceIdxOut) const
605 {
606 return transmissibilities_.thermalHalfTrans(globalSpaceIdxIn,globalSpaceIdxOut);
607 }
608
612 template <class Context>
613 Scalar thermalHalfTransmissibilityIn(const Context& context,
614 unsigned faceIdx,
615 unsigned timeIdx) const
616 {
617 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
618 unsigned toDofLocalIdx = face.exteriorIndex();
619 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransIn;
620 }
621
625 template <class Context>
627 unsigned faceIdx,
628 unsigned timeIdx) const
629 {
630 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
631 unsigned toDofLocalIdx = face.exteriorIndex();
632 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransOut;
633 }
634
638 template <class Context>
640 unsigned boundaryFaceIdx) const
641 {
642 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
643 return transmissibilities_.thermalHalfTransBoundary(elemIdx, boundaryFaceIdx);
644 }
645
649 const typename Vanguard::TransmissibilityType& eclTransmissibilities() const
650 { return transmissibilities_; }
651
652
654 { return tracerModel_; }
655
657 { return tracerModel_; }
658
667 template <class Context>
668 Scalar porosity(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
669 {
670 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
671 return this->porosity(globalSpaceIdx, timeIdx);
672 }
673
680 template <class Context>
681 Scalar dofCenterDepth(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
682 {
683 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
684 return this->dofCenterDepth(globalSpaceIdx);
685 }
686
693 Scalar dofCenterDepth(unsigned globalSpaceIdx) const
694 {
695 return this->simulator().vanguard().cellCenterDepth(globalSpaceIdx);
696 }
697
701 template <class Context>
702 Scalar rockCompressibility(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
703 {
704 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
705 return this->rockCompressibility(globalSpaceIdx);
706 }
707
711 template <class Context>
712 Scalar rockReferencePressure(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
713 {
714 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
715 return this->rockReferencePressure(globalSpaceIdx);
716 }
717
721 template <class Context>
722 const MaterialLawParams& materialLawParams(const Context& context,
723 unsigned spaceIdx, unsigned timeIdx) const
724 {
725 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
726 return this->materialLawParams(globalSpaceIdx);
727 }
728
729 const MaterialLawParams& materialLawParams(unsigned globalDofIdx) const
730 {
731 return materialLawManager_->materialLawParams(globalDofIdx);
732 }
733
734 const MaterialLawParams& materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
735 {
736 return materialLawManager_->materialLawParams(globalDofIdx, facedir);
737 }
738
742 template <class Context>
744 solidEnergyLawParams(const Context& context,
745 unsigned spaceIdx,
746 unsigned timeIdx) const
747 {
748 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
749 return thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
750 }
751
755 template <class Context>
757 thermalConductionLawParams(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
758 {
759 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
760 return thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
761 }
762
769 std::shared_ptr<const EclMaterialLawManager> materialLawManager() const
770 { return materialLawManager_; }
771
772 template <class FluidState>
774 std::array<Evaluation,numPhases> &mobility,
776 FluidState &fluidState,
777 unsigned globalSpaceIdx) const
778 {
779 OPM_TIMEBLOCK_LOCAL(updateRelperms);
780 {
781 // calculate relative permeabilities. note that we store the result into the
782 // mobility_ class attribute. the division by the phase viscosity happens later.
783 const auto& materialParams = materialLawParams(globalSpaceIdx);
784 MaterialLaw::relativePermeabilities(mobility, materialParams, fluidState);
785 Valgrind::CheckDefined(mobility);
786 }
787 if (materialLawManager_->hasDirectionalRelperms()
788 || materialLawManager_->hasDirectionalImbnum())
789 {
790 using Dir = FaceDir::DirEnum;
791 constexpr int ndim = 3;
792 dirMob = std::make_unique<DirectionalMobility<TypeTag, Evaluation>>();
793 Dir facedirs[ndim] = {Dir::XPlus, Dir::YPlus, Dir::ZPlus};
794 for (int i = 0; i<ndim; i++) {
795 const auto& materialParams = materialLawParams(globalSpaceIdx, facedirs[i]);
796 auto& mob_array = dirMob->getArray(i);
797 MaterialLaw::relativePermeabilities(mob_array, materialParams, fluidState);
798 }
799 }
800 }
801
805 std::shared_ptr<EclMaterialLawManager> materialLawManager()
806 { return materialLawManager_; }
807
812 template <class Context>
813 unsigned pvtRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
814 { return pvtRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
815
820 template <class Context>
821 unsigned satnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
822 { return this->satnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
823
828 template <class Context>
829 unsigned miscnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
830 { return this->miscnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
831
836 template <class Context>
837 unsigned plmixnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
838 { return this->plmixnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
839
840 // TODO: polymer related might need to go to the blackoil side
845 template <class Context>
846 Scalar maxPolymerAdsorption(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
847 { return this->maxPolymerAdsorption(context.globalSpaceIndex(spaceIdx, timeIdx)); }
848
852 std::string name() const
853 { return this->simulator().vanguard().caseName(); }
854
858 template <class Context>
859 Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
860 {
861 // use the initial temperature of the DOF if temperature is not a primary
862 // variable
863 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
864 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
865 }
866
867
868 Scalar temperature(unsigned globalDofIdx, unsigned /*timeIdx*/) const
869 {
870 // use the initial temperature of the DOF if temperature is not a primary
871 // variable
872 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
873 }
874
876 solidEnergyLawParams(unsigned globalSpaceIdx,
877 unsigned /*timeIdx*/) const
878 {
879 return this->thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
880 }
882 thermalConductionLawParams(unsigned globalSpaceIdx,
883 unsigned /*timeIdx*/)const
884 {
885 return this->thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
886 }
887
897 Scalar maxOilSaturation(unsigned globalDofIdx) const
898 {
899 if (!this->vapparsActive(this->episodeIndex()))
900 return 0.0;
901
902 return this->maxOilSaturation_[globalDofIdx];
903 }
904
914 void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
915 {
916 if (!this->vapparsActive(this->episodeIndex()))
917 return;
918
919 this->maxOilSaturation_[globalDofIdx] = value;
920 }
921
926 {
927 // Calculate all intensive quantities.
928 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/0);
929
930 // We also need the intensive quantities for timeIdx == 1
931 // corresponding to the start of the current timestep, if we
932 // do not use the storage cache, or if we cannot recycle the
933 // first iteration storage.
934 if (!this->model().enableStorageCache() || !this->recycleFirstIterationStorage()) {
935 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/1);
936 }
937
938 // initialize the wells. Note that this needs to be done after initializing the
939 // intrinsic permeabilities and the after applying the initial solution because
940 // the well model uses these...
941 wellModel_.init();
942
943 aquiferModel_.initialSolutionApplied();
944
945 const bool invalidateFromHyst = updateHysteresis_();
946 if (invalidateFromHyst) {
947 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
948 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
949 }
950 }
951
957 template <class Context>
958 void source(RateVector& rate,
959 const Context& context,
960 unsigned spaceIdx,
961 unsigned timeIdx) const
962 {
963 const unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
964 source(rate, globalDofIdx, timeIdx);
965 }
966
967 void source(RateVector& rate,
968 unsigned globalDofIdx,
969 unsigned timeIdx) const
970 {
971 OPM_TIMEBLOCK_LOCAL(eclProblemSource);
972 rate = 0.0;
973
974 // Add well contribution to source here.
975 wellModel_.computeTotalRatesForDof(rate, globalDofIdx);
976
977 // convert the source term from the total mass rate of the
978 // cell to the one per unit of volume as used by the model.
979 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
980 rate[eqIdx] /= this->model().dofTotalVolume(globalDofIdx);
981
982 Valgrind::CheckDefined(rate[eqIdx]);
983 assert(isfinite(rate[eqIdx]));
984 }
985
986 // Add non-well sources.
987 addToSourceDense(rate, globalDofIdx, timeIdx);
988 }
989
990 virtual void addToSourceDense(RateVector& rate,
991 unsigned globalDofIdx,
992 unsigned timeIdx) const = 0;
993
999 const WellModel& wellModel() const
1000 { return wellModel_; }
1001
1003 { return wellModel_; }
1004
1006 { return aquiferModel_; }
1007
1009 { return aquiferModel_; }
1010
1013
1021 {
1022 OPM_TIMEBLOCK(nexTimeStepSize);
1023 // allow external code to do the timestepping
1024 if (this->nextTimeStepSize_ > 0.0)
1025 return this->nextTimeStepSize_;
1026
1027 const auto& simulator = this->simulator();
1028 int episodeIdx = simulator.episodeIndex();
1029
1030 // for the initial episode, we use a fixed time step size
1031 if (episodeIdx < 0)
1032 return this->initialTimeStepSize_;
1033
1034 // ask the newton method for a suggestion. This suggestion will be based on how
1035 // well the previous time step converged. After that, apply the runtime time
1036 // stepping constraints.
1037 const auto& newtonMethod = this->model().newtonMethod();
1038 return limitNextTimeStepSize_(newtonMethod.suggestTimeStepSize(simulator.timeStepSize()));
1039 }
1040
1046 template <class LhsEval>
1047 LhsEval rockCompPoroMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1048 {
1049 OPM_TIMEBLOCK_LOCAL(rockCompPoroMultiplier);
1050 if (this->rockCompPoroMult_.empty() && this->rockCompPoroMultWc_.empty())
1051 return 1.0;
1052
1053 unsigned tableIdx = 0;
1054 if (!this->rockTableIdx_.empty())
1055 tableIdx = this->rockTableIdx_[elementIdx];
1056
1057 const auto& fs = intQuants.fluidState();
1058 LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
1059 if (!this->minRefPressure_.empty())
1060 // The pore space change is irreversible
1061 effectivePressure =
1062 min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
1063 this->minRefPressure_[elementIdx]);
1064
1065 if (!this->overburdenPressure_.empty())
1066 effectivePressure -= this->overburdenPressure_[elementIdx];
1067
1068
1069 if (!this->rockCompPoroMult_.empty()) {
1070 return this->rockCompPoroMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1071 }
1072
1073 // water compaction
1074 assert(!this->rockCompPoroMultWc_.empty());
1075 LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1076 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1077
1078 return this->rockCompPoroMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1079 }
1080
1086 template <class LhsEval>
1087 LhsEval rockCompTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1088 {
1089 const bool implicit = !this->explicitRockCompaction_;
1090 return implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx)
1091 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1092 }
1093
1094
1098 template <class LhsEval>
1099 LhsEval wellTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1100 {
1101 OPM_TIMEBLOCK_LOCAL(wellTransMultiplier);
1102
1103 const bool implicit = !this->explicitRockCompaction_;
1104 double trans_mult = implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<double>(intQuants, elementIdx)
1105 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1106 trans_mult *= this->simulator().problem().template permFactTransMultiplier<double>(intQuants);
1107
1108 return trans_mult;
1109 }
1110
1111 std::pair<BCType, RateVector> boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
1112 {
1113 OPM_TIMEBLOCK_LOCAL(boundaryCondition);
1115 return { BCType::NONE, RateVector(0.0) };
1116 }
1117 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
1118 const auto& schedule = this->simulator().vanguard().schedule();
1119 if (bcindex_(dir)[globalSpaceIdx] == 0) {
1120 return { BCType::NONE, RateVector(0.0) };
1121 }
1122 if (schedule[this->episodeIndex()].bcprop.size() == 0) {
1123 return { BCType::NONE, RateVector(0.0) };
1124 }
1125 const auto& bc = schedule[this->episodeIndex()].bcprop[bcindex_(dir)[globalSpaceIdx]];
1126 if (bc.bctype!=BCType::RATE) {
1127 return { bc.bctype, RateVector(0.0) };
1128 }
1129
1130 RateVector rate = 0.0;
1131 switch (bc.component) {
1132 case BCComponent::OIL:
1133 rate[Indices::canonicalToActiveComponentIndex(oilCompIdx)] = bc.rate;
1134 break;
1135 case BCComponent::GAS:
1136 rate[Indices::canonicalToActiveComponentIndex(gasCompIdx)] = bc.rate;
1137 break;
1138 case BCComponent::WATER:
1139 rate[Indices::canonicalToActiveComponentIndex(waterCompIdx)] = bc.rate;
1140 break;
1141 case BCComponent::SOLVENT:
1142 this->handleSolventBC(bc, rate);
1143 break;
1144 case BCComponent::POLYMER:
1145 this->handlePolymerBC(bc, rate);
1146 break;
1147 case BCComponent::NONE:
1148 throw std::logic_error("you need to specify the component when RATE type is set in BC");
1149 break;
1150 }
1151 //TODO add support for enthalpy rate
1152 return {bc.bctype, rate};
1153 }
1154
1155
1156 template<class Serializer>
1157 void serializeOp(Serializer& serializer)
1158 {
1159 serializer(static_cast<BaseType&>(*this));
1160 serializer(drift_);
1161 serializer(wellModel_);
1162 serializer(aquiferModel_);
1163 serializer(tracerModel_);
1164 serializer(*materialLawManager_);
1165 }
1166
1167private:
1168 Implementation& asImp_()
1169 { return *static_cast<Implementation *>(this); }
1170
1171 const Implementation& asImp_() const
1172 { return *static_cast<const Implementation *>(this); }
1173
1174protected:
1175 template<class UpdateFunc>
1176 void updateProperty_(const std::string& failureMsg,
1177 UpdateFunc func)
1178 {
1179 OPM_TIMEBLOCK(updateProperty);
1180 const auto& model = this->simulator().model();
1181 const auto& primaryVars = model.solution(/*timeIdx*/0);
1182 const auto& vanguard = this->simulator().vanguard();
1183 std::size_t numGridDof = primaryVars.size();
1185#ifdef _OPENMP
1186#pragma omp parallel for
1187#endif
1188 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
1189 const auto& iq = *model.cachedIntensiveQuantities(dofIdx, /*timeIdx=*/ 0);
1190 func(dofIdx, iq);
1191 }
1192 OPM_END_PARALLEL_TRY_CATCH(failureMsg, vanguard.grid().comm());
1193 }
1194
1196 {
1197 OPM_TIMEBLOCK(updateMaxOilSaturation);
1198 int episodeIdx = this->episodeIndex();
1199
1200 // we use VAPPARS
1201 if (this->vapparsActive(episodeIdx)) {
1202 this->updateProperty_("FlowProblem::updateMaxOilSaturation_() failed:",
1203 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1204 {
1205 this->updateMaxOilSaturation_(compressedDofIdx,iq);
1206 });
1207 return true;
1208 }
1209
1210 return false;
1211 }
1212
1213 bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1214 {
1215 OPM_TIMEBLOCK_LOCAL(updateMaxOilSaturation);
1216 const auto& fs = iq.fluidState();
1217 const Scalar So = decay<Scalar>(fs.saturation(refPressurePhaseIdx_()));
1218 auto& mos = this->maxOilSaturation_;
1219 if(mos[compressedDofIdx] < So){
1220 mos[compressedDofIdx] = So;
1221 return true;
1222 }else{
1223 return false;
1224 }
1225 }
1226
1228 {
1229 OPM_TIMEBLOCK(updateMaxWaterSaturation);
1230 // water compaction is activated in ROCKCOMP
1231 if (this->maxWaterSaturation_.empty())
1232 return false;
1233
1234 this->maxWaterSaturation_[/*timeIdx=*/1] = this->maxWaterSaturation_[/*timeIdx=*/0];
1235 this->updateProperty_("FlowProblem::updateMaxWaterSaturation_() failed:",
1236 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1237 {
1238 this->updateMaxWaterSaturation_(compressedDofIdx,iq);
1239 });
1240 return true;
1241 }
1242
1243
1244 bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1245 {
1246 OPM_TIMEBLOCK_LOCAL(updateMaxWaterSaturation);
1247 const auto& fs = iq.fluidState();
1248 const Scalar Sw = decay<Scalar>(fs.saturation(waterPhaseIdx));
1249 auto& mow = this->maxWaterSaturation_;
1250 if(mow[compressedDofIdx]< Sw){
1251 mow[compressedDofIdx] = Sw;
1252 return true;
1253 }else{
1254 return false;
1255 }
1256 }
1257
1259 {
1260 OPM_TIMEBLOCK(updateMinPressure);
1261 // IRREVERS option is used in ROCKCOMP
1262 if (this->minRefPressure_.empty())
1263 return false;
1264
1265 this->updateProperty_("FlowProblem::updateMinPressure_() failed:",
1266 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1267 {
1268 this->updateMinPressure_(compressedDofIdx,iq);
1269 });
1270 return true;
1271 }
1272
1273 bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities& iq){
1274 OPM_TIMEBLOCK_LOCAL(updateMinPressure);
1275 const auto& fs = iq.fluidState();
1276 const Scalar min_pressure = getValue(fs.pressure(refPressurePhaseIdx_()));
1277 auto& min_pressures = this->minRefPressure_;
1278 if(min_pressures[compressedDofIdx]> min_pressure){
1279 min_pressures[compressedDofIdx] = min_pressure;
1280 return true;
1281 }else{
1282 return false;
1283 }
1284 }
1285
1286 // \brief Function to assign field properties of type double, on the leaf grid view.
1287 //
1288 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1289 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1290 std::function<std::vector<double>(const FieldPropsManager&, const std::string&)>
1292 {
1293 const auto& lookup = this->lookUpData_;
1294 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString)
1295 {
1296 return lookup.assignFieldPropsDoubleOnLeaf(fieldPropManager, propString);
1297 };
1298 }
1299
1300 // \brief Function to assign field properties of type int, unsigned int, ..., on the leaf grid view.
1301 //
1302 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1303 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1304 template<typename IntType>
1305 std::function<std::vector<IntType>(const FieldPropsManager&, const std::string&, bool)>
1307 {
1308 const auto& lookup = this->lookUpData_;
1309 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString, bool needsTranslation)
1310 {
1311 return lookup.template assignFieldPropsIntOnLeaf<IntType>(fieldPropManager, propString, needsTranslation);
1312 };
1313 }
1314
1316 {
1317 OPM_TIMEBLOCK(readMaterialParameters);
1318 const auto& simulator = this->simulator();
1319 const auto& vanguard = simulator.vanguard();
1320 const auto& eclState = vanguard.eclState();
1321
1322 // the PVT and saturation region numbers
1324 this->updatePvtnum_();
1325 this->updateSatnum_();
1326
1327 // the MISC region numbers (solvent model)
1328 this->updateMiscnum_();
1329 // the PLMIX region numbers (polymer model)
1330 this->updatePlmixnum_();
1331
1332 OPM_END_PARALLEL_TRY_CATCH("Invalid region numbers: ", vanguard.gridView().comm());
1334 // porosity
1336 this->referencePorosity_[1] = this->referencePorosity_[0];
1338
1340 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1341 materialLawManager_ = std::make_shared<EclMaterialLawManager>();
1342 materialLawManager_->initFromState(eclState);
1343 materialLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1344 this-> template fieldPropIntTypeOnLeafAssigner_<int>(),
1347 }
1348
1350 {
1351 if constexpr (enableEnergy)
1352 {
1353 const auto& simulator = this->simulator();
1354 const auto& vanguard = simulator.vanguard();
1355 const auto& eclState = vanguard.eclState();
1356
1357 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1358 thermalLawManager_ = std::make_shared<EclThermalLawManager>();
1359 thermalLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1361 this-> template fieldPropIntTypeOnLeafAssigner_<unsigned int>());
1362 }
1363 }
1364
1366 {
1367 const auto& simulator = this->simulator();
1368 const auto& vanguard = simulator.vanguard();
1369 const auto& eclState = vanguard.eclState();
1370
1371 std::size_t numDof = this->model().numGridDof();
1372
1373 this->referencePorosity_[/*timeIdx=*/0].resize(numDof);
1374
1375 const auto& fp = eclState.fieldProps();
1376 const std::vector<double> porvData = this -> fieldPropDoubleOnLeafAssigner_()(fp, "PORV");
1377 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1378 Scalar poreVolume = porvData[dofIdx];
1379
1380 // we define the porosity as the accumulated pore volume divided by the
1381 // geometric volume of the element. Note that -- in pathetic cases -- it can
1382 // be larger than 1.0!
1383 Scalar dofVolume = simulator.model().dofTotalVolume(dofIdx);
1384 assert(dofVolume > 0.0);
1385 this->referencePorosity_[/*timeIdx=*/0][dofIdx] = poreVolume/dofVolume;
1386 }
1387 }
1388
1390 {
1391 // TODO: whether we should move this to FlowProblemBlackoil
1392 const auto& simulator = this->simulator();
1393 const auto& vanguard = simulator.vanguard();
1394 const auto& eclState = vanguard.eclState();
1395
1396 if (eclState.getInitConfig().hasEquil())
1398 else
1400
1401 //initialize min/max values
1402 std::size_t numElems = this->model().numGridDof();
1403 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1404 const auto& fs = asImp_().initialFluidStates()[elemIdx];
1405 if (!this->maxWaterSaturation_.empty() && waterPhaseIdx > -1)
1406 this->maxWaterSaturation_[elemIdx] = std::max(this->maxWaterSaturation_[elemIdx], fs.saturation(waterPhaseIdx));
1407 if (!this->maxOilSaturation_.empty() && oilPhaseIdx > -1)
1408 this->maxOilSaturation_[elemIdx] = std::max(this->maxOilSaturation_[elemIdx], fs.saturation(oilPhaseIdx));
1409 if (!this->minRefPressure_.empty() && refPressurePhaseIdx_() > -1)
1410 this->minRefPressure_[elemIdx] = std::min(this->minRefPressure_[elemIdx], fs.pressure(refPressurePhaseIdx_()));
1411 }
1412 }
1413
1414 virtual void readEquilInitialCondition_() = 0;
1416
1417 // update the hysteresis parameters of the material laws for the whole grid
1419 {
1420 if (!materialLawManager_->enableHysteresis())
1421 return false;
1422
1423 // we need to update the hysteresis data for _all_ elements (i.e., not just the
1424 // interior ones) to avoid desynchronization of the processes in the parallel case!
1425 this->updateProperty_("FlowProblem::updateHysteresis_() failed:",
1426 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1427 {
1428 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1429 });
1430 return true;
1431 }
1432
1433
1434 bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1435 {
1436 OPM_TIMEBLOCK_LOCAL(updateHysteresis_);
1437 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1438 //TODO change materials to give a bool
1439 return true;
1440 }
1441
1442 Scalar getRockCompTransMultVal(std::size_t dofIdx) const
1443 {
1444 if (this->rockCompTransMultVal_.empty())
1445 return 1.0;
1446
1447 return this->rockCompTransMultVal_[dofIdx];
1448 }
1449
1450protected:
1452 {
1453 ConditionalStorage<enableEnergy, Scalar> thermalHalfTransIn;
1454 ConditionalStorage<enableEnergy, Scalar> thermalHalfTransOut;
1455 ConditionalStorage<enableDiffusion, Scalar> diffusivity;
1456 ConditionalStorage<enableDispersion, Scalar> dispersivity;
1458 };
1459
1460 // update the prefetch friendly data object
1462 {
1463 const auto& distFn =
1464 [this](PffDofData_& dofData,
1465 const Stencil& stencil,
1466 unsigned localDofIdx)
1467 -> void
1468 {
1469 const auto& elementMapper = this->model().elementMapper();
1470
1471 unsigned globalElemIdx = elementMapper.index(stencil.entity(localDofIdx));
1472 if (localDofIdx != 0) {
1473 unsigned globalCenterElemIdx = elementMapper.index(stencil.entity(/*dofIdx=*/0));
1474 dofData.transmissibility = transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
1475
1476 if constexpr (enableEnergy) {
1477 *dofData.thermalHalfTransIn = transmissibilities_.thermalHalfTrans(globalCenterElemIdx, globalElemIdx);
1478 *dofData.thermalHalfTransOut = transmissibilities_.thermalHalfTrans(globalElemIdx, globalCenterElemIdx);
1479 }
1480 if constexpr (enableDiffusion)
1481 *dofData.diffusivity = transmissibilities_.diffusivity(globalCenterElemIdx, globalElemIdx);
1482 if (enableDispersion)
1483 dofData.dispersivity = transmissibilities_.dispersivity(globalCenterElemIdx, globalElemIdx);
1484 }
1485 };
1486
1487 pffDofData_.update(distFn);
1488 }
1489
1490 virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart) = 0;
1491
1493 {
1494 const auto& simulator = this->simulator();
1495 const auto& vanguard = simulator.vanguard();
1496 const auto& bcconfig = vanguard.eclState().getSimulationConfig().bcconfig();
1497 if (bcconfig.size() > 0) {
1499
1500 std::size_t numCartDof = vanguard.cartesianSize();
1501 unsigned numElems = vanguard.gridView().size(/*codim=*/0);
1502 std::vector<int> cartesianToCompressedElemIdx(numCartDof, -1);
1503
1504 for (unsigned elemIdx = 0; elemIdx < numElems; ++elemIdx)
1505 cartesianToCompressedElemIdx[vanguard.cartesianIndex(elemIdx)] = elemIdx;
1506
1507 bcindex_.resize(numElems, 0);
1508 auto loopAndApply = [&cartesianToCompressedElemIdx,
1509 &vanguard](const auto& bcface,
1510 auto apply)
1511 {
1512 for (int i = bcface.i1; i <= bcface.i2; ++i) {
1513 for (int j = bcface.j1; j <= bcface.j2; ++j) {
1514 for (int k = bcface.k1; k <= bcface.k2; ++k) {
1515 std::array<int, 3> tmp = {i,j,k};
1516 auto elemIdx = cartesianToCompressedElemIdx[vanguard.cartesianIndex(tmp)];
1517 if (elemIdx >= 0)
1518 apply(elemIdx);
1519 }
1520 }
1521 }
1522 };
1523 for (const auto& bcface : bcconfig) {
1524 std::vector<int>& data = bcindex_(bcface.dir);
1525 const int index = bcface.index;
1526 loopAndApply(bcface,
1527 [&data,index](int elemIdx)
1528 { data[elemIdx] = index; });
1529 }
1530 }
1531 }
1532
1533 // this method applies the runtime constraints specified via the deck and/or command
1534 // line parameters for the size of the next time step.
1536 {
1537 if constexpr (enableExperiments) {
1538 const auto& simulator = this->simulator();
1539 const auto& schedule = simulator.vanguard().schedule();
1540 int episodeIdx = simulator.episodeIndex();
1541
1542 // first thing in the morning, limit the time step size to the maximum size
1543 Scalar maxTimeStepSize = Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60;
1544 int reportStepIdx = std::max(episodeIdx, 0);
1545 if (this->enableTuning_) {
1546 const auto& tuning = schedule[reportStepIdx].tuning();
1547 maxTimeStepSize = tuning.TSMAXZ;
1548 }
1549
1550 dtNext = std::min(dtNext, maxTimeStepSize);
1551
1552 Scalar remainingEpisodeTime =
1553 simulator.episodeStartTime() + simulator.episodeLength()
1554 - (simulator.startTime() + simulator.time());
1555 assert(remainingEpisodeTime >= 0.0);
1556
1557 // if we would have a small amount of time left over in the current episode, make
1558 // two equal time steps instead of a big and a small one
1559 if (remainingEpisodeTime/2.0 < dtNext && dtNext < remainingEpisodeTime*(1.0 - 1e-5))
1560 // note: limiting to the maximum time step size here is probably not strictly
1561 // necessary, but it should not hurt and is more fool-proof
1562 dtNext = std::min(maxTimeStepSize, remainingEpisodeTime/2.0);
1563
1564 if (simulator.episodeStarts()) {
1565 // if a well event occurred, respect the limit for the maximum time step after
1566 // that, too
1567 const auto& events = simulator.vanguard().schedule()[reportStepIdx].events();
1568 bool wellEventOccured =
1569 events.hasEvent(ScheduleEvents::NEW_WELL)
1570 || events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE)
1571 || events.hasEvent(ScheduleEvents::INJECTION_UPDATE)
1572 || events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
1573 if (episodeIdx >= 0 && wellEventOccured && this->maxTimeStepAfterWellEvent_ > 0)
1574 dtNext = std::min(dtNext, this->maxTimeStepAfterWellEvent_);
1575 }
1576 }
1577
1578 return dtNext;
1579 }
1580
1582 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1583 return oilPhaseIdx;
1584 }
1585 else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
1586 return gasPhaseIdx;
1587 }
1588 else {
1589 return waterPhaseIdx;
1590 }
1591 }
1592
1594 {
1595 const auto& model = this->simulator().model();
1596 std::size_t numGridDof = this->model().numGridDof();
1597 this->rockCompTransMultVal_.resize(numGridDof, 1.0);
1598 for (std::size_t elementIdx = 0; elementIdx < numGridDof; ++elementIdx) {
1599 const auto& iq = *model.cachedIntensiveQuantities(elementIdx, /*timeIdx=*/ 0);
1600 Scalar trans_mult = computeRockCompTransMultiplier_<Scalar>(iq, elementIdx);
1601 this->rockCompTransMultVal_[elementIdx] = trans_mult;
1602 }
1603 }
1604
1610 template <class LhsEval>
1611 LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1612 {
1613 OPM_TIMEBLOCK_LOCAL(computeRockCompTransMultiplier);
1614 if (this->rockCompTransMult_.empty() && this->rockCompTransMultWc_.empty())
1615 return 1.0;
1616
1617 unsigned tableIdx = 0;
1618 if (!this->rockTableIdx_.empty())
1619 tableIdx = this->rockTableIdx_[elementIdx];
1620
1621 const auto& fs = intQuants.fluidState();
1622 LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
1623
1624 if (!this->minRefPressure_.empty())
1625 // The pore space change is irreversible
1626 effectivePressure =
1627 min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
1628 this->minRefPressure_[elementIdx]);
1629
1630 if (!this->overburdenPressure_.empty())
1631 effectivePressure -= this->overburdenPressure_[elementIdx];
1632
1633 if (!this->rockCompTransMult_.empty())
1634 return this->rockCompTransMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1635
1636 // water compaction
1637 assert(!this->rockCompTransMultWc_.empty());
1638 LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1639 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1640
1641 return this->rockCompTransMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1642 }
1643
1644 typename Vanguard::TransmissibilityType transmissibilities_;
1645
1646 std::shared_ptr<EclMaterialLawManager> materialLawManager_;
1647 std::shared_ptr<EclThermalLawManager> thermalLawManager_;
1648
1651
1654
1656
1657
1660
1661 template<class T>
1662 struct BCData
1663 {
1664 std::array<std::vector<T>,6> data;
1665
1666 void resize(std::size_t size, T defVal)
1667 {
1668 for (auto& d : data)
1669 d.resize(size, defVal);
1670 }
1671
1672 const std::vector<T>& operator()(FaceDir::DirEnum dir) const
1673 {
1674 if (dir == FaceDir::DirEnum::Unknown)
1675 throw std::runtime_error("Tried to access BC data for the 'Unknown' direction");
1676 int idx = 0;
1677 int div = static_cast<int>(dir);
1678 while ((div /= 2) >= 1)
1679 ++idx;
1680 assert(idx >= 0 && idx <= 5);
1681 return data[idx];
1682 }
1683
1684 std::vector<T>& operator()(FaceDir::DirEnum dir)
1685 {
1686 return const_cast<std::vector<T>&>(std::as_const(*this)(dir));
1687 }
1688 };
1689
1690 virtual void handleSolventBC(const BCProp::BCFace&, RateVector&) const = 0;
1691
1692 virtual void handlePolymerBC(const BCProp::BCFace&, RateVector&) const = 0;
1693
1697 bool first_step_ = true;
1698
1699};
1700
1701} // namespace Opm
1702
1703#endif // OPM_FLOW_PROBLEM_HPP
#define OPM_END_PARALLEL_TRY_CATCH(prefix, comm)
Catch exception and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:192
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowGenericProblem.hpp:70
Scalar maxPolymerAdsorption(unsigned elemIdx) const
Returns the max polymer adsorption value.
Definition: FlowGenericProblem_impl.hpp:743
unsigned pvtRegionIndex(unsigned elemIdx) const
Returns the index the relevant PVT region given a cell index.
Definition: FlowGenericProblem_impl.hpp:702
std::function< unsigned(unsigned)> lookupIdxOnLevelZeroAssigner_()
Definition: FlowGenericProblem.hpp:383
Scalar rockReferencePressure(unsigned globalSpaceIdx) const
Definition: FlowGenericProblem_impl.hpp:324
std::vector< TabulatedTwoDFunction > rockCompPoroMultWc_
Definition: FlowGenericProblem.hpp:351
Scalar porosity(unsigned globalSpaceIdx, unsigned timeIdx) const
Direct indexed access to the porosity.
Definition: FlowGenericProblem_impl.hpp:339
Scalar rockCompressibility(unsigned globalSpaceIdx) const
Definition: FlowGenericProblem_impl.hpp:309
unsigned miscnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant MISC region given a cell index.
Definition: FlowGenericProblem_impl.hpp:722
unsigned satnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant saturation function region given a cell index.
Definition: FlowGenericProblem_impl.hpp:712
void beginTimeStep_(bool enableExperiments, int episodeIdx, int timeStepIndex, Scalar startTime, Scalar time, Scalar timeStepSize, Scalar endTime)
Definition: FlowGenericProblem_impl.hpp:458
std::vector< TabulatedTwoDFunction > rockCompTransMultWc_
Definition: FlowGenericProblem.hpp:352
unsigned plmixnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant PLMIXNUM (for polymer module) region given a cell index.
Definition: FlowGenericProblem_impl.hpp:732
std::array< std::vector< Scalar >, 2 > referencePorosity_
Definition: FlowGenericProblem.hpp:342
bool beginEpisode_(bool enableExperiments, int episodeIdx)
Definition: FlowGenericProblem_impl.hpp:421
bool shouldWriteOutput() const
Always returns true. The ecl output writer takes care of the rest.
Definition: FlowGenericProblem.hpp:300
static std::string helpPreamble(int, const char **argv)
Definition: FlowGenericProblem_impl.hpp:116
bool shouldWriteRestartFile() const
Returns true if an eWoms restart file should be written to disk.
Definition: FlowGenericProblem.hpp:309
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblem.hpp:91
const WellModel & wellModel() const
Returns a reference to the ECL well manager used by the problem.
Definition: FlowProblem.hpp:999
GetPropType< TypeTag, Properties::Evaluation > Evaluation
Definition: FlowProblem.hpp:153
Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
Direct access to the transmissibility between two elements.
Definition: FlowProblem.hpp:536
static int handlePositionalParameter(std::function< void(const std::string &, const std::string &)> addKey, std::set< std::string > &seenParams, std::string &errorMsg, int, const char **argv, int paramIdx, int)
Handles positional command line parameters.
Definition: FlowProblem.hpp:189
@ numComponents
Definition: FlowProblem.hpp:113
bool nonTrivialBoundaryConditions() const
Definition: FlowProblem.hpp:1011
std::function< std::vector< IntType >(const FieldPropsManager &, const std::string &, bool)> fieldPropIntTypeOnLeafAssigner_()
Definition: FlowProblem.hpp:1306
const DimMatrix & intrinsicPermeability(unsigned globalElemIdx) const
This method returns the intrinsic permeability tensor given a global element index.
Definition: FlowProblem.hpp:518
@ enableSolvent
Definition: FlowProblem.hpp:127
@ enablePolymer
Definition: FlowProblem.hpp:124
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:813
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:668
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition: FlowProblem.hpp:395
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:104
GetPropType< TypeTag, Properties::DofMapper > DofMapper
Definition: FlowProblem.hpp:152
@ enableMICP
Definition: FlowProblem.hpp:123
@ waterCompIdx
Definition: FlowProblem.hpp:139
@ enableDiffusion
Definition: FlowProblem.hpp:117
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: FlowProblem.hpp:98
typename EclThermalLawManager::SolidEnergyLawParams SolidEnergyLawParams
Definition: FlowProblem.hpp:149
bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1434
GetPropType< TypeTag, Properties::BaseProblem > ParentType
Definition: FlowProblem.hpp:95
GetPropType< TypeTag, Properties::EqVector > EqVector
Definition: FlowProblem.hpp:103
unsigned satnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:821
virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart)=0
@ numEq
Definition: FlowProblem.hpp:111
Scalar thermalHalfTransmissibilityOut(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:626
@ enableEnergy
Definition: FlowProblem.hpp:119
bool first_step_
Definition: FlowProblem.hpp:1697
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:145
const ThermalConductionLawParams & thermalConductionLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:882
AquiferModel aquiferModel_
Definition: FlowProblem.hpp:1653
GlobalEqVector drift_
Definition: FlowProblem.hpp:1650
bool updateMinPressure_()
Definition: FlowProblem.hpp:1258
std::function< std::vector< double >(const FieldPropsManager &, const std::string &)> fieldPropDoubleOnLeafAssigner_()
Definition: FlowProblem.hpp:1291
bool enableDriftCompensation_
Definition: FlowProblem.hpp:1649
Scalar transmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:583
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:142
void updateReferencePorosity_()
Definition: FlowProblem.hpp:1365
Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn, const unsigned globalSpaceIdxOut) const
Definition: FlowProblem.hpp:603
BCData< int > bcindex_
Definition: FlowProblem.hpp:1694
GetPropType< TypeTag, Properties::TracerModel > TracerModel
Definition: FlowProblem.hpp:163
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:129
@ dim
Definition: FlowProblem.hpp:107
Scalar rockCompressibility(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:702
bool updateMaxWaterSaturation_()
Definition: FlowProblem.hpp:1227
Utility::CopyablePtr< DirectionalMobility< TypeTag, Evaluation > > DirectionalMobilityPtr
Definition: FlowProblem.hpp:164
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: FlowProblem.hpp:160
Scalar maxOilSaturation(unsigned globalDofIdx) const
Returns an element's historic maximum oil phase saturation that was observed during the simulation.
Definition: FlowProblem.hpp:897
bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1273
std::string name() const
The problem name.
Definition: FlowProblem.hpp:852
int episodeIndex() const
Definition: FlowProblem.hpp:297
LhsEval rockCompTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1087
void endIteration()
Called by the simulator after each Newton-Raphson iteration.
Definition: FlowProblem.hpp:405
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:154
FlowProblem(Simulator &simulator)
Definition: FlowProblem.hpp:208
@ enableTemperature
Definition: FlowProblem.hpp:128
@ gasPhaseIdx
Definition: FlowProblem.hpp:131
@ oilCompIdx
Definition: FlowProblem.hpp:138
virtual void handleSolventBC(const BCProp::BCFace &, RateVector &) const =0
@ enableExperiments
Definition: FlowProblem.hpp:120
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: FlowProblem.hpp:102
@ dimWorld
Definition: FlowProblem.hpp:108
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: FlowProblem.hpp:143
const Vanguard::TransmissibilityType & eclTransmissibilities() const
Return a reference to the object that handles the "raw" transmissibilities.
Definition: FlowProblem.hpp:649
void source(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:967
virtual void readEquilInitialCondition_()=0
Scalar nextTimeStepSize() const
Propose the size of the next time step to the simulator.
Definition: FlowProblem.hpp:1020
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1611
std::pair< BCType, RateVector > boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
Definition: FlowProblem.hpp:1111
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:126
typename EclMaterialLawManager::MaterialLawParams MaterialLawParams
Definition: FlowProblem.hpp:148
virtual ~FlowProblem()=default
PffGridVector< GridView, Stencil, PffDofData_, DofMapper > pffDofData_
Definition: FlowProblem.hpp:1658
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:773
@ enablePolymerMolarWeight
Definition: FlowProblem.hpp:125
const ThermalConductionLawParams & thermalConductionLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:757
Scalar transmissibility(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:525
@ enableConvectiveMixing
Definition: FlowProblem.hpp:115
TracerModel tracerModel_
Definition: FlowProblem.hpp:1659
virtual void handlePolymerBC(const BCProp::BCFace &, RateVector &) const =0
const TracerModel & tracerModel() const
Definition: FlowProblem.hpp:653
WellModel wellModel_
Definition: FlowProblem.hpp:1652
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:305
const SolidEnergyLawParams & solidEnergyLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:876
@ numPhases
Definition: FlowProblem.hpp:112
Scalar getRockCompTransMultVal(std::size_t dofIdx) const
Definition: FlowProblem.hpp:1442
Scalar dofCenterDepth(unsigned globalSpaceIdx) const
Direct indexed acces to the depth of an degree of freedom [m].
Definition: FlowProblem.hpp:693
typename GetProp< TypeTag, Properties::MaterialLaw >::EclMaterialLawManager EclMaterialLawManager
Definition: FlowProblem.hpp:146
LhsEval wellTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Return the well transmissibility multiplier due to rock changues.
Definition: FlowProblem.hpp:1099
bool explicitRockCompaction_
Definition: FlowProblem.hpp:1696
const MaterialLawParams & materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
Definition: FlowProblem.hpp:734
std::shared_ptr< const EclMaterialLawManager > materialLawManager() const
Returns the ECL material law manager.
Definition: FlowProblem.hpp:769
unsigned plmixnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:837
const MaterialLawParams & materialLawParams(unsigned globalDofIdx) const
Definition: FlowProblem.hpp:729
Scalar temperature(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:859
Scalar thermalHalfTransmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:639
unsigned miscnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:829
void updateRockCompTransMultVal_()
Definition: FlowProblem.hpp:1593
std::shared_ptr< EclThermalLawManager > thermalLawManager_
Definition: FlowProblem.hpp:1647
Scalar limitNextTimeStepSize_(Scalar dtNext) const
Definition: FlowProblem.hpp:1535
GetPropType< TypeTag, Properties::Stencil > Stencil
Definition: FlowProblem.hpp:100
typename GetProp< TypeTag, Properties::SolidEnergyLaw >::EclThermalLawManager EclThermalLawManager
Definition: FlowProblem.hpp:147
virtual void readExplicitInitialCondition_()=0
GetPropType< TypeTag, Properties::WellModel > WellModel
Definition: FlowProblem.hpp:156
const DimMatrix & intrinsicPermeability(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:504
bool updateHysteresis_()
Definition: FlowProblem.hpp:1418
void readThermalParameters_()
Definition: FlowProblem.hpp:1349
@ enableFoam
Definition: FlowProblem.hpp:122
void serializeOp(Serializer &serializer)
Definition: FlowProblem.hpp:1157
@ enableBrine
Definition: FlowProblem.hpp:116
Scalar maxPolymerAdsorption(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the max polymer adsorption value.
Definition: FlowProblem.hpp:846
void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:362
Scalar thermalHalfTransmissibilityIn(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:613
void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
Sets an element's maximum oil phase saturation observed during the simulation.
Definition: FlowProblem.hpp:914
@ waterPhaseIdx
Definition: FlowProblem.hpp:133
const SolidEnergyLawParams & solidEnergyLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Return the parameters for the energy storage law of the rock.
Definition: FlowProblem.hpp:744
typename EclThermalLawManager::ThermalConductionLawParams ThermalConductionLawParams
Definition: FlowProblem.hpp:150
Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:563
AquiferModel & mutableAquiferModel()
Definition: FlowProblem.hpp:1008
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:155
Scalar temperature(unsigned globalDofIdx, unsigned) const
Definition: FlowProblem.hpp:868
GetPropType< TypeTag, Properties::AquiferModel > AquiferModel
Definition: FlowProblem.hpp:157
bool enableVtkOutput_
Definition: FlowProblem.hpp:1655
std::shared_ptr< EclMaterialLawManager > materialLawManager_
Definition: FlowProblem.hpp:1646
WellModel & wellModel()
Definition: FlowProblem.hpp:1002
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowProblem.hpp:99
void updateProperty_(const std::string &failureMsg, UpdateFunc func)
Definition: FlowProblem.hpp:1176
bool updateMaxOilSaturation_()
Definition: FlowProblem.hpp:1195
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblem.hpp:178
void updatePffDofData_()
Definition: FlowProblem.hpp:1461
virtual void endEpisode()
Called by the simulator after the end of an episode.
Definition: FlowProblem.hpp:466
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:141
bool nonTrivialBoundaryConditions_
Definition: FlowProblem.hpp:1695
@ enableExtbo
Definition: FlowProblem.hpp:121
@ gasCompIdx
Definition: FlowProblem.hpp:137
GetPropType< TypeTag, Properties::Problem > Implementation
Definition: FlowProblem.hpp:96
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1492
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1644
void source(RateVector &rate, const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Evaluate the source term for all phases within a given sub-control-volume.
Definition: FlowProblem.hpp:958
virtual void addToSourceDense(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const =0
virtual void endTimeStep()
Called by the simulator after each time integration.
Definition: FlowProblem.hpp:415
virtual void readInitialCondition_()
Definition: FlowProblem.hpp:1389
virtual void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblem.hpp:925
const MaterialLawParams & materialLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:722
Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:570
void writeOutput(bool verbose=true)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:489
@ oilPhaseIdx
Definition: FlowProblem.hpp:132
Scalar diffusivity(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:545
void deserialize(Restarter &res)
This method restores the complete state of the problem and its sub-objects from disk.
Definition: FlowProblem.hpp:271
std::shared_ptr< EclMaterialLawManager > materialLawManager()
Definition: FlowProblem.hpp:805
@ enableDispersion
Definition: FlowProblem.hpp:118
bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1244
typename GridView::template Codim< 0 >::Entity Element
Definition: FlowProblem.hpp:144
Scalar dofCenterDepth(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the depth of an degree of freedom [m].
Definition: FlowProblem.hpp:681
bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1213
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: FlowProblem.hpp:101
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: FlowProblem.hpp:151
const AquiferModel & aquiferModel() const
Definition: FlowProblem.hpp:1005
int refPressurePhaseIdx_() const
Definition: FlowProblem.hpp:1581
LhsEval rockCompPoroMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the porosity multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1047
Scalar transmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:593
TracerModel & tracerModel()
Definition: FlowProblem.hpp:656
void readMaterialParameters_()
Definition: FlowProblem.hpp:1315
void serialize(Restarter &res)
This method writes the complete state of the problem and its subobjects to disk.
Definition: FlowProblem.hpp:290
MathToolbox< Evaluation > Toolbox
Definition: FlowProblem.hpp:159
Scalar rockReferencePressure(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:712
Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:556
void prefetch(const Element &elem) const
Definition: FlowProblem.hpp:256
A random-access container which stores data attached to a grid's degrees of freedom in a prefetch fri...
Definition: pffgridvector.hh:49
Definition: RelpermDiagnostics.hpp:50
void diagnosis(const EclipseState &eclState, const CartesianIndexMapper &cartesianIndexMapper)
This file contains definitions related to directional mobilities.
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilboundaryratevector.hh:37
int eclPositionalParameter(std::function< void(const std::string &, const std::string &)> addKey, std::set< std::string > &seenParams, std::string &errorMsg, const char **argv, int paramIdx)
Definition: FlowGenericProblem_impl.hpp:48
void eclBroadcast(Parallel::Communication, T &)
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:235
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:226
Definition: FlowProblem.hpp:1663
const std::vector< T > & operator()(FaceDir::DirEnum dir) const
Definition: FlowProblem.hpp:1672
void resize(std::size_t size, T defVal)
Definition: FlowProblem.hpp:1666
std::vector< T > & operator()(FaceDir::DirEnum dir)
Definition: FlowProblem.hpp:1684
std::array< std::vector< T >, 6 > data
Definition: FlowProblem.hpp:1664
Definition: FlowProblem.hpp:1452
ConditionalStorage< enableEnergy, Scalar > thermalHalfTransIn
Definition: FlowProblem.hpp:1453
ConditionalStorage< enableEnergy, Scalar > thermalHalfTransOut
Definition: FlowProblem.hpp:1454
ConditionalStorage< enableDiffusion, Scalar > diffusivity
Definition: FlowProblem.hpp:1455
ConditionalStorage< enableDispersion, Scalar > dispersivity
Definition: FlowProblem.hpp:1456
Scalar transmissibility
Definition: FlowProblem.hpp:1457