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/Schedule/Schedule.hpp>
41#include <opm/input/eclipse/Units/Units.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
67
71
72#include <opm/utility/CopyablePtr.hpp>
73
74#include <algorithm>
75#include <cstddef>
76#include <functional>
77#include <set>
78#include <stdexcept>
79#include <string>
80#include <vector>
81
82namespace Opm {
83
90template <class TypeTag>
91class FlowProblem : public GetPropType<TypeTag, Properties::BaseProblem>
92 , public FlowGenericProblem<GetPropType<TypeTag, Properties::GridView>,
93 GetPropType<TypeTag, Properties::FluidSystem>>
94{
95protected:
100
108
109 // Grid and world dimension
110 enum { dim = GridView::dimension };
111 enum { dimWorld = GridView::dimensionworld };
112
113 // copy some indices for convenience
114 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
115 enum { numPhases = FluidSystem::numPhases };
116 enum { numComponents = FluidSystem::numComponents };
117
118 enum { enableConvectiveMixing = getPropValue<TypeTag, Properties::EnableConvectiveMixing>() };
119 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
120 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
121 enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
122 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
123 enum { enableExperiments = getPropValue<TypeTag, Properties::EnableExperiments>() };
124 enum { enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>() };
125 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
126 enum { enableMICP = getPropValue<TypeTag, Properties::EnableMICP>() };
127 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
128 enum { enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>() };
129 enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
130 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
131 enum { enableTemperature = getPropValue<TypeTag, Properties::EnableTemperature>() };
132 enum { enableThermalFluxBoundaries = getPropValue<TypeTag, Properties::EnableThermalFluxBoundaries>() };
133
134 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
135 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
136 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
137
138 // TODO: later, gasCompIdx, oilCompIdx and waterCompIdx should go to the FlowProblemBlackoil in the future
139 // we do not want them in the compositional setting
140 enum { gasCompIdx = FluidSystem::gasCompIdx };
141 enum { oilCompIdx = FluidSystem::oilCompIdx };
142 enum { waterCompIdx = FluidSystem::waterCompIdx };
143
147 using Element = typename GridView::template Codim<0>::Entity;
151 using MaterialLawParams = typename EclMaterialLawManager::MaterialLawParams;
152 using SolidEnergyLawParams = typename EclThermalLawManager::SolidEnergyLawParams;
153 using ThermalConductionLawParams = typename EclThermalLawManager::ThermalConductionLawParams;
161
162 using Toolbox = MathToolbox<Evaluation>;
163 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
164
165
167 using DirectionalMobilityPtr = Utility::CopyablePtr<DirectionalMobility<TypeTag>>;
168
169public:
175 using BaseType::porosity;
176
180 static void registerParameters()
181 {
182 ParentType::registerParameters();
183
184 registerFlowProblemParameters<Scalar>();
185 }
186
190 static int handlePositionalParameter(std::function<void(const std::string&,
191 const std::string&)> addKey,
192 std::set<std::string>& seenParams,
193 std::string& errorMsg,
194 int,
195 const char** argv,
196 int paramIdx,
197 int)
198 {
199 return detail::eclPositionalParameter(addKey,
200 seenParams,
201 errorMsg,
202 argv,
203 paramIdx);
204 }
205
209 explicit FlowProblem(Simulator& simulator)
210 : ParentType(simulator)
211 , BaseType(simulator.vanguard().eclState(),
212 simulator.vanguard().schedule(),
213 simulator.vanguard().gridView())
214 , transmissibilities_(simulator.vanguard().eclState(),
215 simulator.vanguard().gridView(),
216 simulator.vanguard().cartesianIndexMapper(),
217 simulator.vanguard().grid(),
218 simulator.vanguard().cellCentroids(),
222 , wellModel_(simulator)
223 , aquiferModel_(simulator)
224 , pffDofData_(simulator.gridView(), this->elementMapper())
225 , tracerModel_(simulator)
226 {
227 if (! Parameters::Get<Parameters::CheckSatfuncConsistency>()) {
228 // User did not enable the "new" saturation function consistency
229 // check module. Run the original checker instead. This is a
230 // temporary measure.
231 RelpermDiagnostics relpermDiagnostics{};
232 relpermDiagnostics.diagnosis(simulator.vanguard().eclState(),
233 simulator.vanguard().levelCartesianIndexMapper());
234 }
235 }
236
237 virtual ~FlowProblem() = default;
238
239 void prefetch(const Element& elem) const
240 { this->pffDofData_.prefetch(elem); }
241
253 template <class Restarter>
254 void deserialize(Restarter& res)
255 {
256 // reload the current episode/report step from the deck
257 this->beginEpisode();
258
259 // deserialize the wells
260 wellModel_.deserialize(res);
261
262 // deserialize the aquifer
263 aquiferModel_.deserialize(res);
264 }
265
272 template <class Restarter>
273 void serialize(Restarter& res)
274 {
275 wellModel_.serialize(res);
276
277 aquiferModel_.serialize(res);
278 }
279
280 int episodeIndex() const
281 {
282 return std::max(this->simulator().episodeIndex(), 0);
283 }
284
288 virtual void beginEpisode()
289 {
290 OPM_TIMEBLOCK(beginEpisode);
291 // Proceed to the next report step
292 auto& simulator = this->simulator();
293 int episodeIdx = simulator.episodeIndex();
294 auto& eclState = simulator.vanguard().eclState();
295 const auto& schedule = simulator.vanguard().schedule();
296 const auto& events = schedule[episodeIdx].events();
297
298 if (episodeIdx >= 0 && events.hasEvent(ScheduleEvents::GEO_MODIFIER)) {
299 // bring the contents of the keywords to the current state of the SCHEDULE
300 // section.
301 //
302 // TODO (?): make grid topology changes possible (depending on what exactly
303 // has changed, the grid may need be re-created which has some serious
304 // implications on e.g., the solution of the simulation.)
305 const auto& miniDeck = schedule[episodeIdx].geo_keywords();
306 const auto& cc = simulator.vanguard().grid().comm();
307 eclState.apply_schedule_keywords( miniDeck );
308 eclBroadcast(cc, eclState.getTransMult() );
309
310 // Re-ordering in case of ALUGrid
311 std::function<unsigned int(unsigned int)> equilGridToGrid = [&simulator](unsigned int i) {
312 return simulator.vanguard().gridEquilIdxToGridIdx(i);
313 };
314
315 // re-compute all quantities which may possibly be affected.
316 using TransUpdateQuantities = typename Vanguard::TransmissibilityType::TransUpdateQuantities;
317 transmissibilities_.update(true, TransUpdateQuantities::All, equilGridToGrid);
318 this->referencePorosity_[1] = this->referencePorosity_[0];
321 this->model().linearizer().updateDiscretizationParameters();
322 }
323
324 bool tuningEvent = this->beginEpisode_(enableExperiments, this->episodeIndex());
325
326 // set up the wells for the next episode.
327 wellModel_.beginEpisode();
328
329 // set up the aquifers for the next episode.
330 aquiferModel_.beginEpisode();
331
332 // set the size of the initial time step of the episode
333 Scalar dt = limitNextTimeStepSize_(simulator.episodeLength());
334 // negative value of initialTimeStepSize_ indicates no active limit from TSINIT or NEXTSTEP
335 if ( (episodeIdx == 0 || tuningEvent) && this->initialTimeStepSize_ > 0)
336 // allow the size of the initial time step to be set via an external parameter
337 // if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
338 dt = std::min(dt, this->initialTimeStepSize_);
339 simulator.setTimeStepSize(dt);
340 }
341
346 {
347 OPM_TIMEBLOCK(beginTimeStep);
348 const int episodeIdx = this->episodeIndex();
349 const int timeStepSize = this->simulator().timeStepSize();
350
352 episodeIdx,
353 this->simulator().timeStepIndex(),
354 this->simulator().startTime(),
355 this->simulator().time(),
356 timeStepSize,
357 this->simulator().endTime());
358
359 // update maximum water saturation and minimum pressure
360 // used when ROCKCOMP is activated
361 // Do not update max RS first step after a restart
362 this->updateExplicitQuantities_(episodeIdx, timeStepSize, first_step_ && (episodeIdx > 0));
363 first_step_ = false;
364
366 this->model().linearizer().updateBoundaryConditionData();
367 }
368
369 wellModel_.beginTimeStep();
370 aquiferModel_.beginTimeStep();
371 tracerModel_.beginTimeStep();
372
373 }
374
379 {
380 OPM_TIMEBLOCK(beginIteration);
381 wellModel_.beginIteration();
382 aquiferModel_.beginIteration();
383 }
384
389 {
390 OPM_TIMEBLOCK(endIteration);
391 wellModel_.endIteration();
392 aquiferModel_.endIteration();
393 }
394
398 virtual void endTimeStep()
399 {
400 OPM_TIMEBLOCK(endTimeStep);
401
402#ifndef NDEBUG
403 if constexpr (getPropValue<TypeTag, Properties::EnableDebuggingChecks>()) {
404 // in debug mode, we don't care about performance, so we check
405 // if the model does the right thing (i.e., the mass change
406 // inside the whole reservoir must be equivalent to the fluxes
407 // over the grid's boundaries plus the source rates specified by
408 // the problem).
409 const int rank = this->simulator().gridView().comm().rank();
410 if (rank == 0) {
411 std::cout << "checking conservativeness of solution\n";
412 }
413
414 this->model().checkConservativeness(/*tolerance=*/-1, /*verbose=*/true);
415 if (rank == 0) {
416 std::cout << "solution is sufficiently conservative\n";
417 }
418 }
419#endif // NDEBUG
420
421 auto& simulator = this->simulator();
422 simulator.setTimeStepIndex(simulator.timeStepIndex()+1);
423
424 this->wellModel_.endTimeStep();
425 this->aquiferModel_.endTimeStep();
426 this->tracerModel_.endTimeStep();
427
428 // Compute flux for output
429 this->model().linearizer().updateFlowsInfo();
430
431 if (this->enableDriftCompensation_) {
432 OPM_TIMEBLOCK(driftCompansation);
433
434 const auto& residual = this->model().linearizer().residual();
435
436 for (unsigned globalDofIdx = 0; globalDofIdx < residual.size(); globalDofIdx ++) {
437 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(globalDofIdx);
438 this->drift_[sfcdofIdx] = residual[sfcdofIdx] * simulator.timeStepSize();
439
440 if constexpr (getPropValue<TypeTag, Properties::UseVolumetricResidual>()) {
441 this->drift_[sfcdofIdx] *= this->model().dofTotalVolume(sfcdofIdx);
442 }
443 }
444 }
445 }
446
450 virtual void endEpisode()
451 {
452 const int episodeIdx = this->episodeIndex();
453
454 this->wellModel_.endEpisode();
455 this->aquiferModel_.endEpisode();
456
457 const auto& schedule = this->simulator().vanguard().schedule();
458
459 // End simulation when completed.
460 if (episodeIdx + 1 >= static_cast<int>(schedule.size()) - 1) {
461 this->simulator().setFinished(true);
462 return;
463 }
464
465 // Otherwise, start next episode (report step).
466 this->simulator().startNextEpisode(schedule.stepLength(episodeIdx + 1));
467 }
468
473 virtual void writeOutput(bool verbose)
474 {
475 OPM_TIMEBLOCK(problemWriteOutput);
476 // use the generic code to prepare the output fields and to
477 // write the desired VTK files.
478 if (Parameters::Get<Parameters::EnableWriteAllSolutions>() ||
479 this->simulator().episodeWillBeOver()) {
480 ParentType::writeOutput(verbose);
481 }
482 }
483
487 template <class Context>
488 const DimMatrix& intrinsicPermeability(const Context& context,
489 unsigned spaceIdx,
490 unsigned timeIdx) const
491 {
492 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
493 return transmissibilities_.permeability(globalSpaceIdx);
494 }
495
502 const DimMatrix& intrinsicPermeability(unsigned globalElemIdx) const
503 { return transmissibilities_.permeability(globalElemIdx); }
504
508 template <class Context>
509 Scalar transmissibility(const Context& context,
510 [[maybe_unused]] unsigned fromDofLocalIdx,
511 unsigned toDofLocalIdx) const
512 {
513 assert(fromDofLocalIdx == 0);
514 return pffDofData_.get(context.element(), toDofLocalIdx).transmissibility;
515 }
516
520 Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
521 {
522 return transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
523 }
524
528 template <class Context>
529 Scalar diffusivity(const Context& context,
530 [[maybe_unused]] unsigned fromDofLocalIdx,
531 unsigned toDofLocalIdx) const
532 {
533 assert(fromDofLocalIdx == 0);
534 return *pffDofData_.get(context.element(), toDofLocalIdx).diffusivity;
535 }
536
540 Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
541 return transmissibilities_.diffusivity(globalCellIn, globalCellOut);
542 }
543
547 Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
548 return transmissibilities_.dispersivity(globalCellIn, globalCellOut);
549 }
550
554 Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx,
555 const unsigned boundaryFaceIdx) const
556 {
557 return transmissibilities_.thermalTransmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
558 }
559
560
561
562
566 template <class Context>
567 Scalar transmissibilityBoundary(const Context& elemCtx,
568 unsigned boundaryFaceIdx) const
569 {
570 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
571 return transmissibilities_.transmissibilityBoundary(elemIdx, boundaryFaceIdx);
572 }
573
577 Scalar transmissibilityBoundary(const unsigned globalSpaceIdx,
578 const unsigned boundaryFaceIdx) const
579 {
580 return transmissibilities_.transmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
581 }
582
583
587 Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn,
588 const unsigned globalSpaceIdxOut) const
589 {
590 return transmissibilities_.thermalHalfTrans(globalSpaceIdxIn,globalSpaceIdxOut);
591 }
592
596 template <class Context>
597 Scalar thermalHalfTransmissibilityIn(const Context& context,
598 unsigned faceIdx,
599 unsigned timeIdx) const
600 {
601 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
602 unsigned toDofLocalIdx = face.exteriorIndex();
603 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransIn;
604 }
605
609 template <class Context>
611 unsigned faceIdx,
612 unsigned timeIdx) const
613 {
614 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
615 unsigned toDofLocalIdx = face.exteriorIndex();
616 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransOut;
617 }
618
622 template <class Context>
624 unsigned boundaryFaceIdx) const
625 {
626 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
627 return transmissibilities_.thermalHalfTransBoundary(elemIdx, boundaryFaceIdx);
628 }
629
633 const typename Vanguard::TransmissibilityType& eclTransmissibilities() const
634 { return transmissibilities_; }
635
636
638 { return tracerModel_; }
639
641 { return tracerModel_; }
642
651 template <class Context>
652 Scalar porosity(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
653 {
654 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
655 return this->porosity(globalSpaceIdx, timeIdx);
656 }
657
664 template <class Context>
665 Scalar dofCenterDepth(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
666 {
667 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
668 return this->dofCenterDepth(globalSpaceIdx);
669 }
670
677 Scalar dofCenterDepth(unsigned globalSpaceIdx) const
678 {
679 return this->simulator().vanguard().cellCenterDepth(globalSpaceIdx);
680 }
681
685 template <class Context>
686 Scalar rockCompressibility(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
687 {
688 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
689 return this->rockCompressibility(globalSpaceIdx);
690 }
691
695 template <class Context>
696 Scalar rockReferencePressure(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
697 {
698 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
699 return rockReferencePressure(globalSpaceIdx);
700 }
701
705 Scalar rockReferencePressure(unsigned globalSpaceIdx) const
706 {
707 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
708 if (rock_config.store()) {
709 return asImp_().initialFluidState(globalSpaceIdx).pressure(refPressurePhaseIdx_());
710 }
711 else {
712 if (this->rockParams_.empty())
713 return 1e5;
714
715 unsigned tableIdx = 0;
716 if (!this->rockTableIdx_.empty()) {
717 tableIdx = this->rockTableIdx_[globalSpaceIdx];
718 }
719 return this->rockParams_[tableIdx].referencePressure;
720 }
721 }
722
726 template <class Context>
727 const MaterialLawParams& materialLawParams(const Context& context,
728 unsigned spaceIdx, unsigned timeIdx) const
729 {
730 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
731 return this->materialLawParams(globalSpaceIdx);
732 }
733
734 const MaterialLawParams& materialLawParams(unsigned globalDofIdx) const
735 {
736 return materialLawManager_->materialLawParams(globalDofIdx);
737 }
738
739 const MaterialLawParams& materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
740 {
741 return materialLawManager_->materialLawParams(globalDofIdx, facedir);
742 }
743
747 template <class Context>
749 solidEnergyLawParams(const Context& context,
750 unsigned spaceIdx,
751 unsigned timeIdx) const
752 {
753 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
754 return thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
755 }
756
760 template <class Context>
762 thermalConductionLawParams(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
763 {
764 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
765 return thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
766 }
767
774 std::shared_ptr<const EclMaterialLawManager> materialLawManager() const
775 { return materialLawManager_; }
776
777 template <class FluidState>
779 std::array<Evaluation,numPhases> &mobility,
781 FluidState &fluidState,
782 unsigned globalSpaceIdx) const
783 {
784 OPM_TIMEBLOCK_LOCAL(updateRelperms);
785 {
786 // calculate relative permeabilities. note that we store the result into the
787 // mobility_ class attribute. the division by the phase viscosity happens later.
788 const auto& materialParams = materialLawParams(globalSpaceIdx);
789 MaterialLaw::relativePermeabilities(mobility, materialParams, fluidState);
790 Valgrind::CheckDefined(mobility);
791 }
792 if (materialLawManager_->hasDirectionalRelperms()
793 || materialLawManager_->hasDirectionalImbnum())
794 {
795 using Dir = FaceDir::DirEnum;
796 constexpr int ndim = 3;
797 dirMob = std::make_unique<DirectionalMobility<TypeTag>>();
798 Dir facedirs[ndim] = {Dir::XPlus, Dir::YPlus, Dir::ZPlus};
799 for (int i = 0; i<ndim; i++) {
800 const auto& materialParams = materialLawParams(globalSpaceIdx, facedirs[i]);
801 auto& mob_array = dirMob->getArray(i);
802 MaterialLaw::relativePermeabilities(mob_array, materialParams, fluidState);
803 }
804 }
805 }
806
810 std::shared_ptr<EclMaterialLawManager> materialLawManager()
811 { return materialLawManager_; }
812
817 template <class Context>
818 unsigned pvtRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
819 { return pvtRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
820
825 template <class Context>
826 unsigned satnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
827 { return this->satnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
828
833 template <class Context>
834 unsigned miscnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
835 { return this->miscnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
836
841 template <class Context>
842 unsigned plmixnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
843 { return this->plmixnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
844
845 // TODO: polymer related might need to go to the blackoil side
850 template <class Context>
851 Scalar maxPolymerAdsorption(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
852 { return this->maxPolymerAdsorption(context.globalSpaceIndex(spaceIdx, timeIdx)); }
853
857 std::string name() const
858 { return this->simulator().vanguard().caseName(); }
859
863 template <class Context>
864 Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
865 {
866 // use the initial temperature of the DOF if temperature is not a primary
867 // variable
868 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
869 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
870 }
871
872
873 Scalar temperature(unsigned globalDofIdx, unsigned /*timeIdx*/) const
874 {
875 // use the initial temperature of the DOF if temperature is not a primary
876 // variable
877 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
878 }
879
881 solidEnergyLawParams(unsigned globalSpaceIdx,
882 unsigned /*timeIdx*/) const
883 {
884 return this->thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
885 }
887 thermalConductionLawParams(unsigned globalSpaceIdx,
888 unsigned /*timeIdx*/)const
889 {
890 return this->thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
891 }
892
902 Scalar maxOilSaturation(unsigned globalDofIdx) const
903 {
904 if (!this->vapparsActive(this->episodeIndex()))
905 return 0.0;
906
907 return this->maxOilSaturation_[globalDofIdx];
908 }
909
919 void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
920 {
921 if (!this->vapparsActive(this->episodeIndex()))
922 return;
923
924 this->maxOilSaturation_[globalDofIdx] = value;
925 }
926
931 {
932 // Calculate all intensive quantities.
933 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/0);
934
935 // We also need the intensive quantities for timeIdx == 1
936 // corresponding to the start of the current timestep, if we
937 // do not use the storage cache, or if we cannot recycle the
938 // first iteration storage.
939 if (!this->model().enableStorageCache() || !this->recycleFirstIterationStorage()) {
940 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/1);
941 }
942
943 // initialize the wells. Note that this needs to be done after initializing the
944 // intrinsic permeabilities and the after applying the initial solution because
945 // the well model uses these...
946 wellModel_.init();
947
948 aquiferModel_.initialSolutionApplied();
949
950 const bool invalidateFromHyst = updateHysteresis_();
951 if (invalidateFromHyst) {
952 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
953 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
954 }
955 }
956
962 template <class Context>
963 void source(RateVector& rate,
964 const Context& context,
965 unsigned spaceIdx,
966 unsigned timeIdx) const
967 {
968 const unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
969 source(rate, globalDofIdx, timeIdx);
970 }
971
972 void source(RateVector& rate,
973 unsigned globalDofIdx,
974 unsigned timeIdx) const
975 {
976 OPM_TIMEBLOCK_LOCAL(eclProblemSource);
977 rate = 0.0;
978
979 // Add well contribution to source here.
980 wellModel_.computeTotalRatesForDof(rate, globalDofIdx);
981
982 // convert the source term from the total mass rate of the
983 // cell to the one per unit of volume as used by the model.
984 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
985 rate[eqIdx] /= this->model().dofTotalVolume(globalDofIdx);
986
987 Valgrind::CheckDefined(rate[eqIdx]);
988 assert(isfinite(rate[eqIdx]));
989 }
990
991 // Add non-well sources.
992 addToSourceDense(rate, globalDofIdx, timeIdx);
993 }
994
995 virtual void addToSourceDense(RateVector& rate,
996 unsigned globalDofIdx,
997 unsigned timeIdx) const = 0;
998
1004 const WellModel& wellModel() const
1005 { return wellModel_; }
1006
1008 { return wellModel_; }
1009
1011 { return aquiferModel_; }
1012
1014 { return aquiferModel_; }
1015
1018
1026 {
1027 OPM_TIMEBLOCK(nexTimeStepSize);
1028 // allow external code to do the timestepping
1029 if (this->nextTimeStepSize_ > 0.0)
1030 return this->nextTimeStepSize_;
1031
1032 const auto& simulator = this->simulator();
1033 int episodeIdx = simulator.episodeIndex();
1034
1035 // for the initial episode, we use a fixed time step size
1036 if (episodeIdx < 0)
1037 return this->initialTimeStepSize_;
1038
1039 // ask the newton method for a suggestion. This suggestion will be based on how
1040 // well the previous time step converged. After that, apply the runtime time
1041 // stepping constraints.
1042 const auto& newtonMethod = this->model().newtonMethod();
1043 return limitNextTimeStepSize_(newtonMethod.suggestTimeStepSize(simulator.timeStepSize()));
1044 }
1045
1051 template <class LhsEval>
1052 LhsEval rockCompPoroMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1053 {
1054 OPM_TIMEBLOCK_LOCAL(rockCompPoroMultiplier);
1055 if (this->rockCompPoroMult_.empty() && this->rockCompPoroMultWc_.empty())
1056 return 1.0;
1057
1058 unsigned tableIdx = 0;
1059 if (!this->rockTableIdx_.empty())
1060 tableIdx = this->rockTableIdx_[elementIdx];
1061
1062 const auto& fs = intQuants.fluidState();
1063 LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
1064 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1065 if (!this->minRefPressure_.empty())
1066 // The pore space change is irreversible
1067 effectivePressure =
1068 min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
1069 this->minRefPressure_[elementIdx]);
1070
1071 if (!this->overburdenPressure_.empty())
1072 effectivePressure -= this->overburdenPressure_[elementIdx];
1073
1074 if (rock_config.store()) {
1075 effectivePressure -= asImp_().initialFluidState(elementIdx).pressure(refPressurePhaseIdx_());
1076 }
1077
1078 if (!this->rockCompPoroMult_.empty()) {
1079 return this->rockCompPoroMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1080 }
1081
1082 // water compaction
1083 assert(!this->rockCompPoroMultWc_.empty());
1084 LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1085 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1086
1087 return this->rockCompPoroMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1088 }
1089
1095 template <class LhsEval>
1096 LhsEval rockCompTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1097 {
1098 const bool implicit = !this->explicitRockCompaction_;
1099 return implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx)
1100 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1101 }
1102
1103
1107 template <class LhsEval>
1108 LhsEval wellTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1109 {
1110 OPM_TIMEBLOCK_LOCAL(wellTransMultiplier);
1111
1112 const bool implicit = !this->explicitRockCompaction_;
1113 double trans_mult = implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<double>(intQuants, elementIdx)
1114 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1115 trans_mult *= this->simulator().problem().template permFactTransMultiplier<double>(intQuants, elementIdx);
1116
1117 return trans_mult;
1118 }
1119
1120 std::pair<BCType, RateVector> boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
1121 {
1122 OPM_TIMEBLOCK_LOCAL(boundaryCondition);
1124 return { BCType::NONE, RateVector(0.0) };
1125 }
1126 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
1127 const auto& schedule = this->simulator().vanguard().schedule();
1128 if (bcindex_(dir)[globalSpaceIdx] == 0) {
1129 return { BCType::NONE, RateVector(0.0) };
1130 }
1131 if (schedule[this->episodeIndex()].bcprop.size() == 0) {
1132 return { BCType::NONE, RateVector(0.0) };
1133 }
1134 const auto& bc = schedule[this->episodeIndex()].bcprop[bcindex_(dir)[globalSpaceIdx]];
1135 if (bc.bctype!=BCType::RATE) {
1136 return { bc.bctype, RateVector(0.0) };
1137 }
1138
1139 RateVector rate = 0.0;
1140 switch (bc.component) {
1141 case BCComponent::OIL:
1142 rate[Indices::canonicalToActiveComponentIndex(oilCompIdx)] = bc.rate;
1143 break;
1144 case BCComponent::GAS:
1145 rate[Indices::canonicalToActiveComponentIndex(gasCompIdx)] = bc.rate;
1146 break;
1147 case BCComponent::WATER:
1148 rate[Indices::canonicalToActiveComponentIndex(waterCompIdx)] = bc.rate;
1149 break;
1150 case BCComponent::SOLVENT:
1151 this->handleSolventBC(bc, rate);
1152 break;
1153 case BCComponent::POLYMER:
1154 this->handlePolymerBC(bc, rate);
1155 break;
1156 case BCComponent::MICR:
1157 this->handleMicrBC(bc, rate);
1158 break;
1159 case BCComponent::OXYG:
1160 this->handleOxygBC(bc, rate);
1161 break;
1162 case BCComponent::UREA:
1163 this->handleUreaBC(bc, rate);
1164 break;
1165 case BCComponent::NONE:
1166 throw std::logic_error("you need to specify the component when RATE type is set in BC");
1167 break;
1168 }
1169 //TODO add support for enthalpy rate
1170 return {bc.bctype, rate};
1171 }
1172
1173
1174 template<class Serializer>
1175 void serializeOp(Serializer& serializer)
1176 {
1177 serializer(static_cast<BaseType&>(*this));
1178 serializer(drift_);
1179 serializer(wellModel_);
1180 serializer(aquiferModel_);
1181 serializer(tracerModel_);
1182 serializer(*materialLawManager_);
1183 }
1184
1185private:
1186 Implementation& asImp_()
1187 { return *static_cast<Implementation *>(this); }
1188
1189 const Implementation& asImp_() const
1190 { return *static_cast<const Implementation *>(this); }
1191
1192protected:
1193 template<class UpdateFunc>
1194 void updateProperty_(const std::string& failureMsg,
1195 UpdateFunc func)
1196 {
1197 OPM_TIMEBLOCK(updateProperty);
1198 const auto& model = this->simulator().model();
1199 const auto& primaryVars = model.solution(/*timeIdx*/0);
1200 const auto& vanguard = this->simulator().vanguard();
1201 std::size_t numGridDof = primaryVars.size();
1203#ifdef _OPENMP
1204#pragma omp parallel for
1205#endif
1206 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
1207 const auto& iq = *model.cachedIntensiveQuantities(dofIdx, /*timeIdx=*/ 0);
1208 func(dofIdx, iq);
1209 }
1210 OPM_END_PARALLEL_TRY_CATCH(failureMsg, vanguard.grid().comm());
1211 }
1212
1214 {
1215 OPM_TIMEBLOCK(updateMaxOilSaturation);
1216 int episodeIdx = this->episodeIndex();
1217
1218 // we use VAPPARS
1219 if (this->vapparsActive(episodeIdx)) {
1220 this->updateProperty_("FlowProblem::updateMaxOilSaturation_() failed:",
1221 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1222 {
1223 this->updateMaxOilSaturation_(compressedDofIdx,iq);
1224 });
1225 return true;
1226 }
1227
1228 return false;
1229 }
1230
1231 bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1232 {
1233 OPM_TIMEBLOCK_LOCAL(updateMaxOilSaturation);
1234 const auto& fs = iq.fluidState();
1235 const Scalar So = decay<Scalar>(fs.saturation(refPressurePhaseIdx_()));
1236 auto& mos = this->maxOilSaturation_;
1237 if(mos[compressedDofIdx] < So){
1238 mos[compressedDofIdx] = So;
1239 return true;
1240 }else{
1241 return false;
1242 }
1243 }
1244
1246 {
1247 OPM_TIMEBLOCK(updateMaxWaterSaturation);
1248 // water compaction is activated in ROCKCOMP
1249 if (this->maxWaterSaturation_.empty())
1250 return false;
1251
1252 this->maxWaterSaturation_[/*timeIdx=*/1] = this->maxWaterSaturation_[/*timeIdx=*/0];
1253 this->updateProperty_("FlowProblem::updateMaxWaterSaturation_() failed:",
1254 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1255 {
1256 this->updateMaxWaterSaturation_(compressedDofIdx,iq);
1257 });
1258 return true;
1259 }
1260
1261
1262 bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1263 {
1264 OPM_TIMEBLOCK_LOCAL(updateMaxWaterSaturation);
1265 const auto& fs = iq.fluidState();
1266 const Scalar Sw = decay<Scalar>(fs.saturation(waterPhaseIdx));
1267 auto& mow = this->maxWaterSaturation_;
1268 if(mow[compressedDofIdx]< Sw){
1269 mow[compressedDofIdx] = Sw;
1270 return true;
1271 }else{
1272 return false;
1273 }
1274 }
1275
1277 {
1278 OPM_TIMEBLOCK(updateMinPressure);
1279 // IRREVERS option is used in ROCKCOMP
1280 if (this->minRefPressure_.empty())
1281 return false;
1282
1283 this->updateProperty_("FlowProblem::updateMinPressure_() failed:",
1284 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1285 {
1286 this->updateMinPressure_(compressedDofIdx,iq);
1287 });
1288 return true;
1289 }
1290
1291 bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities& iq){
1292 OPM_TIMEBLOCK_LOCAL(updateMinPressure);
1293 const auto& fs = iq.fluidState();
1294 const Scalar min_pressure = getValue(fs.pressure(refPressurePhaseIdx_()));
1295 auto& min_pressures = this->minRefPressure_;
1296 if(min_pressures[compressedDofIdx]> min_pressure){
1297 min_pressures[compressedDofIdx] = min_pressure;
1298 return true;
1299 }else{
1300 return false;
1301 }
1302 }
1303
1304 // \brief Function to assign field properties of type double, on the leaf grid view.
1305 //
1306 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1307 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1308 std::function<std::vector<double>(const FieldPropsManager&, const std::string&)>
1310 {
1311 const auto& lookup = this->lookUpData_;
1312 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString)
1313 {
1314 return lookup.assignFieldPropsDoubleOnLeaf(fieldPropManager, propString);
1315 };
1316 }
1317
1318 // \brief Function to assign field properties of type int, unsigned int, ..., on the leaf grid view.
1319 //
1320 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1321 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1322 template<typename IntType>
1323 std::function<std::vector<IntType>(const FieldPropsManager&, const std::string&, bool)>
1325 {
1326 const auto& lookup = this->lookUpData_;
1327 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString, bool needsTranslation)
1328 {
1329 return lookup.template assignFieldPropsIntOnLeaf<IntType>(fieldPropManager, propString, needsTranslation);
1330 };
1331 }
1332
1334 {
1335 OPM_TIMEBLOCK(readMaterialParameters);
1336 const auto& simulator = this->simulator();
1337 const auto& vanguard = simulator.vanguard();
1338 const auto& eclState = vanguard.eclState();
1339
1340 // the PVT and saturation region numbers
1342 this->updatePvtnum_();
1343 this->updateSatnum_();
1344
1345 // the MISC region numbers (solvent model)
1346 this->updateMiscnum_();
1347 // the PLMIX region numbers (polymer model)
1348 this->updatePlmixnum_();
1349
1350 OPM_END_PARALLEL_TRY_CATCH("Invalid region numbers: ", vanguard.gridView().comm());
1352 // porosity
1354 this->referencePorosity_[1] = this->referencePorosity_[0];
1356
1358 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1359 materialLawManager_ = std::make_shared<EclMaterialLawManager>();
1360 materialLawManager_->initFromState(eclState);
1361 materialLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1362 this-> template fieldPropIntTypeOnLeafAssigner_<int>(),
1365 }
1366
1368 {
1369 if constexpr (enableEnergy)
1370 {
1371 const auto& simulator = this->simulator();
1372 const auto& vanguard = simulator.vanguard();
1373 const auto& eclState = vanguard.eclState();
1374
1375 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1376 thermalLawManager_ = std::make_shared<EclThermalLawManager>();
1377 thermalLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1379 this-> template fieldPropIntTypeOnLeafAssigner_<unsigned int>());
1380 }
1381 }
1382
1384 {
1385 const auto& simulator = this->simulator();
1386 const auto& vanguard = simulator.vanguard();
1387 const auto& eclState = vanguard.eclState();
1388
1389 std::size_t numDof = this->model().numGridDof();
1390
1391 this->referencePorosity_[/*timeIdx=*/0].resize(numDof);
1392
1393 const auto& fp = eclState.fieldProps();
1394 const std::vector<double> porvData = this -> fieldPropDoubleOnLeafAssigner_()(fp, "PORV");
1395 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1396 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(dofIdx);
1397 Scalar poreVolume = porvData[dofIdx];
1398
1399 // we define the porosity as the accumulated pore volume divided by the
1400 // geometric volume of the element. Note that -- in pathetic cases -- it can
1401 // be larger than 1.0!
1402 Scalar dofVolume = simulator.model().dofTotalVolume(sfcdofIdx);
1403 assert(dofVolume > 0.0);
1404 this->referencePorosity_[/*timeIdx=*/0][sfcdofIdx] = poreVolume/dofVolume;
1405 }
1406 }
1407
1409 {
1410 // TODO: whether we should move this to FlowProblemBlackoil
1411 const auto& simulator = this->simulator();
1412 const auto& vanguard = simulator.vanguard();
1413 const auto& eclState = vanguard.eclState();
1414
1415 if (eclState.getInitConfig().hasEquil())
1417 else
1419
1420 //initialize min/max values
1421 std::size_t numElems = this->model().numGridDof();
1422 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1423 const auto& fs = asImp_().initialFluidStates()[elemIdx];
1424 if (!this->maxWaterSaturation_.empty() && waterPhaseIdx > -1)
1425 this->maxWaterSaturation_[elemIdx] = std::max(this->maxWaterSaturation_[elemIdx], fs.saturation(waterPhaseIdx));
1426 if (!this->maxOilSaturation_.empty() && oilPhaseIdx > -1)
1427 this->maxOilSaturation_[elemIdx] = std::max(this->maxOilSaturation_[elemIdx], fs.saturation(oilPhaseIdx));
1428 if (!this->minRefPressure_.empty() && refPressurePhaseIdx_() > -1)
1429 this->minRefPressure_[elemIdx] = std::min(this->minRefPressure_[elemIdx], fs.pressure(refPressurePhaseIdx_()));
1430 }
1431 }
1432
1433 virtual void readEquilInitialCondition_() = 0;
1435
1436 // update the hysteresis parameters of the material laws for the whole grid
1438 {
1439 if (!materialLawManager_->enableHysteresis())
1440 return false;
1441
1442 // we need to update the hysteresis data for _all_ elements (i.e., not just the
1443 // interior ones) to avoid desynchronization of the processes in the parallel case!
1444 this->updateProperty_("FlowProblem::updateHysteresis_() failed:",
1445 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1446 {
1447 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1448 });
1449 return true;
1450 }
1451
1452
1453 bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1454 {
1455 OPM_TIMEBLOCK_LOCAL(updateHysteresis_);
1456 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1457 //TODO change materials to give a bool
1458 return true;
1459 }
1460
1461 Scalar getRockCompTransMultVal(std::size_t dofIdx) const
1462 {
1463 if (this->rockCompTransMultVal_.empty())
1464 return 1.0;
1465
1466 return this->rockCompTransMultVal_[dofIdx];
1467 }
1468
1469protected:
1471 {
1472 ConditionalStorage<enableEnergy, Scalar> thermalHalfTransIn;
1473 ConditionalStorage<enableEnergy, Scalar> thermalHalfTransOut;
1474 ConditionalStorage<enableDiffusion, Scalar> diffusivity;
1475 ConditionalStorage<enableDispersion, Scalar> dispersivity;
1477 };
1478
1479 // update the prefetch friendly data object
1481 {
1482 const auto& distFn =
1483 [this](PffDofData_& dofData,
1484 const Stencil& stencil,
1485 unsigned localDofIdx)
1486 -> void
1487 {
1488 const auto& elementMapper = this->model().elementMapper();
1489
1490 unsigned globalElemIdx = elementMapper.index(stencil.entity(localDofIdx));
1491 if (localDofIdx != 0) {
1492 unsigned globalCenterElemIdx = elementMapper.index(stencil.entity(/*dofIdx=*/0));
1493 dofData.transmissibility = transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
1494
1495 if constexpr (enableEnergy) {
1496 *dofData.thermalHalfTransIn = transmissibilities_.thermalHalfTrans(globalCenterElemIdx, globalElemIdx);
1497 *dofData.thermalHalfTransOut = transmissibilities_.thermalHalfTrans(globalElemIdx, globalCenterElemIdx);
1498 }
1499 if constexpr (enableDiffusion)
1500 *dofData.diffusivity = transmissibilities_.diffusivity(globalCenterElemIdx, globalElemIdx);
1501 if (enableDispersion)
1502 dofData.dispersivity = transmissibilities_.dispersivity(globalCenterElemIdx, globalElemIdx);
1503 }
1504 };
1505
1506 pffDofData_.update(distFn);
1507 }
1508
1509 virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart) = 0;
1510
1512 {
1513 const auto& simulator = this->simulator();
1514 const auto& vanguard = simulator.vanguard();
1515 const auto& bcconfig = vanguard.eclState().getSimulationConfig().bcconfig();
1516 if (bcconfig.size() > 0) {
1518
1519 std::size_t numCartDof = vanguard.cartesianSize();
1520 unsigned numElems = vanguard.gridView().size(/*codim=*/0);
1521 std::vector<int> cartesianToCompressedElemIdx(numCartDof, -1);
1522
1523 for (unsigned elemIdx = 0; elemIdx < numElems; ++elemIdx)
1524 cartesianToCompressedElemIdx[vanguard.cartesianIndex(elemIdx)] = elemIdx;
1525
1526 bcindex_.resize(numElems, 0);
1527 auto loopAndApply = [&cartesianToCompressedElemIdx,
1528 &vanguard](const auto& bcface,
1529 auto apply)
1530 {
1531 for (int i = bcface.i1; i <= bcface.i2; ++i) {
1532 for (int j = bcface.j1; j <= bcface.j2; ++j) {
1533 for (int k = bcface.k1; k <= bcface.k2; ++k) {
1534 std::array<int, 3> tmp = {i,j,k};
1535 auto elemIdx = cartesianToCompressedElemIdx[vanguard.cartesianIndex(tmp)];
1536 if (elemIdx >= 0)
1537 apply(elemIdx);
1538 }
1539 }
1540 }
1541 };
1542 for (const auto& bcface : bcconfig) {
1543 std::vector<int>& data = bcindex_(bcface.dir);
1544 const int index = bcface.index;
1545 loopAndApply(bcface,
1546 [&data,index](int elemIdx)
1547 { data[elemIdx] = index; });
1548 }
1549 }
1550 }
1551
1552 // this method applies the runtime constraints specified via the deck and/or command
1553 // line parameters for the size of the next time step.
1555 {
1556 if constexpr (enableExperiments) {
1557 const auto& simulator = this->simulator();
1558 const auto& schedule = simulator.vanguard().schedule();
1559 int episodeIdx = simulator.episodeIndex();
1560
1561 // first thing in the morning, limit the time step size to the maximum size
1562 Scalar maxTimeStepSize = Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60;
1563 int reportStepIdx = std::max(episodeIdx, 0);
1564 if (this->enableTuning_) {
1565 const auto& tuning = schedule[reportStepIdx].tuning();
1566 maxTimeStepSize = tuning.TSMAXZ;
1567 }
1568
1569 dtNext = std::min(dtNext, maxTimeStepSize);
1570
1571 Scalar remainingEpisodeTime =
1572 simulator.episodeStartTime() + simulator.episodeLength()
1573 - (simulator.startTime() + simulator.time());
1574 assert(remainingEpisodeTime >= 0.0);
1575
1576 // if we would have a small amount of time left over in the current episode, make
1577 // two equal time steps instead of a big and a small one
1578 if (remainingEpisodeTime/2.0 < dtNext && dtNext < remainingEpisodeTime*(1.0 - 1e-5))
1579 // note: limiting to the maximum time step size here is probably not strictly
1580 // necessary, but it should not hurt and is more fool-proof
1581 dtNext = std::min(maxTimeStepSize, remainingEpisodeTime/2.0);
1582
1583 if (simulator.episodeStarts()) {
1584 // if a well event occurred, respect the limit for the maximum time step after
1585 // that, too
1586 const auto& events = simulator.vanguard().schedule()[reportStepIdx].events();
1587 bool wellEventOccured =
1588 events.hasEvent(ScheduleEvents::NEW_WELL)
1589 || events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE)
1590 || events.hasEvent(ScheduleEvents::INJECTION_UPDATE)
1591 || events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
1592 if (episodeIdx >= 0 && wellEventOccured && this->maxTimeStepAfterWellEvent_ > 0)
1593 dtNext = std::min(dtNext, this->maxTimeStepAfterWellEvent_);
1594 }
1595 }
1596
1597 return dtNext;
1598 }
1599
1601 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1602 return oilPhaseIdx;
1603 }
1604 else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
1605 return gasPhaseIdx;
1606 }
1607 else {
1608 return waterPhaseIdx;
1609 }
1610 }
1611
1613 {
1614 const auto& model = this->simulator().model();
1615 std::size_t numGridDof = this->model().numGridDof();
1616 this->rockCompTransMultVal_.resize(numGridDof, 1.0);
1617 for (std::size_t elementIdx = 0; elementIdx < numGridDof; ++elementIdx) {
1618 const auto& iq = *model.cachedIntensiveQuantities(elementIdx, /*timeIdx=*/ 0);
1619 Scalar trans_mult = computeRockCompTransMultiplier_<Scalar>(iq, elementIdx);
1620 this->rockCompTransMultVal_[elementIdx] = trans_mult;
1621 }
1622 }
1623
1629 template <class LhsEval>
1630 LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1631 {
1632 OPM_TIMEBLOCK_LOCAL(computeRockCompTransMultiplier);
1633 if (this->rockCompTransMult_.empty() && this->rockCompTransMultWc_.empty())
1634 return 1.0;
1635
1636 unsigned tableIdx = 0;
1637 if (!this->rockTableIdx_.empty())
1638 tableIdx = this->rockTableIdx_[elementIdx];
1639
1640 const auto& fs = intQuants.fluidState();
1641 LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
1642 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1643 if (!this->minRefPressure_.empty())
1644 // The pore space change is irreversible
1645 effectivePressure =
1646 min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
1647 this->minRefPressure_[elementIdx]);
1648
1649 if (!this->overburdenPressure_.empty())
1650 effectivePressure -= this->overburdenPressure_[elementIdx];
1651
1652 if (rock_config.store()) {
1653 effectivePressure -= asImp_().initialFluidState(elementIdx).pressure(refPressurePhaseIdx_());
1654 }
1655
1656 if (!this->rockCompTransMult_.empty())
1657 return this->rockCompTransMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1658
1659 // water compaction
1660 assert(!this->rockCompTransMultWc_.empty());
1661 LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1662 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1663
1664 return this->rockCompTransMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1665 }
1666
1667 typename Vanguard::TransmissibilityType transmissibilities_;
1668
1669 std::shared_ptr<EclMaterialLawManager> materialLawManager_;
1670 std::shared_ptr<EclThermalLawManager> thermalLawManager_;
1671
1673
1676
1679
1680 template<class T>
1681 struct BCData
1682 {
1683 std::array<std::vector<T>,6> data;
1684
1685 void resize(std::size_t size, T defVal)
1686 {
1687 for (auto& d : data)
1688 d.resize(size, defVal);
1689 }
1690
1691 const std::vector<T>& operator()(FaceDir::DirEnum dir) const
1692 {
1693 if (dir == FaceDir::DirEnum::Unknown)
1694 throw std::runtime_error("Tried to access BC data for the 'Unknown' direction");
1695 int idx = 0;
1696 int div = static_cast<int>(dir);
1697 while ((div /= 2) >= 1)
1698 ++idx;
1699 assert(idx >= 0 && idx <= 5);
1700 return data[idx];
1701 }
1702
1703 std::vector<T>& operator()(FaceDir::DirEnum dir)
1704 {
1705 return const_cast<std::vector<T>&>(std::as_const(*this)(dir));
1706 }
1707 };
1708
1709 virtual void handleSolventBC(const BCProp::BCFace&, RateVector&) const = 0;
1710
1711 virtual void handlePolymerBC(const BCProp::BCFace&, RateVector&) const = 0;
1712
1713 virtual void handleMicrBC(const BCProp::BCFace&, RateVector&) const = 0;
1714
1715 virtual void handleOxygBC(const BCProp::BCFace&, RateVector&) const = 0;
1716
1717 virtual void handleUreaBC(const BCProp::BCFace&, RateVector&) const = 0;
1718
1721 bool first_step_ = true;
1722};
1723
1724} // namespace Opm
1725
1726#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:61
Scalar maxPolymerAdsorption(unsigned elemIdx) const
Returns the max polymer adsorption value.
Definition: FlowGenericProblem_impl.hpp:733
unsigned pvtRegionIndex(unsigned elemIdx) const
Returns the index the relevant PVT region given a cell index.
Definition: FlowGenericProblem_impl.hpp:692
std::function< unsigned(unsigned)> lookupIdxOnLevelZeroAssigner_()
Definition: FlowGenericProblem.hpp:363
std::vector< TabulatedTwoDFunction > rockCompPoroMultWc_
Definition: FlowGenericProblem.hpp:328
Scalar porosity(unsigned globalSpaceIdx, unsigned timeIdx) const
Direct indexed access to the porosity.
Definition: FlowGenericProblem_impl.hpp:329
Scalar rockCompressibility(unsigned globalSpaceIdx) const
Definition: FlowGenericProblem_impl.hpp:314
unsigned miscnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant MISC region given a cell index.
Definition: FlowGenericProblem_impl.hpp:712
unsigned satnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant saturation function region given a cell index.
Definition: FlowGenericProblem_impl.hpp:702
void beginTimeStep_(bool enableExperiments, int episodeIdx, int timeStepIndex, Scalar startTime, Scalar time, Scalar timeStepSize, Scalar endTime)
Definition: FlowGenericProblem_impl.hpp:448
std::vector< TabulatedTwoDFunction > rockCompTransMultWc_
Definition: FlowGenericProblem.hpp:329
unsigned plmixnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant PLMIXNUM (for polymer module) region given a cell index.
Definition: FlowGenericProblem_impl.hpp:722
std::array< std::vector< Scalar >, 2 > referencePorosity_
Definition: FlowGenericProblem.hpp:319
bool beginEpisode_(bool enableExperiments, int episodeIdx)
Definition: FlowGenericProblem_impl.hpp:411
bool shouldWriteOutput() const
Always returns true. The ecl output writer takes care of the rest.
Definition: FlowGenericProblem.hpp:277
static std::string helpPreamble(int, const char **argv)
Definition: FlowGenericProblem_impl.hpp:111
bool shouldWriteRestartFile() const
Returns true if an eWoms restart file should be written to disk.
Definition: FlowGenericProblem.hpp:286
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblem.hpp:94
const WellModel & wellModel() const
Returns a reference to the ECL well manager used by the problem.
Definition: FlowProblem.hpp:1004
GetPropType< TypeTag, Properties::Evaluation > Evaluation
Definition: FlowProblem.hpp:156
Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
Direct access to the transmissibility between two elements.
Definition: FlowProblem.hpp:520
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:190
@ dimWorld
Definition: FlowProblem.hpp:111
bool nonTrivialBoundaryConditions() const
Definition: FlowProblem.hpp:1016
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:473
std::function< std::vector< IntType >(const FieldPropsManager &, const std::string &, bool)> fieldPropIntTypeOnLeafAssigner_()
Definition: FlowProblem.hpp:1324
const DimMatrix & intrinsicPermeability(unsigned globalElemIdx) const
This method returns the intrinsic permeability tensor given a global element index.
Definition: FlowProblem.hpp:502
@ oilCompIdx
Definition: FlowProblem.hpp:141
@ enablePolymer
Definition: FlowProblem.hpp:127
@ enableExtbo
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:818
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:652
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition: FlowProblem.hpp:378
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:107
GetPropType< TypeTag, Properties::DofMapper > DofMapper
Definition: FlowProblem.hpp:155
@ enableExperiments
Definition: FlowProblem.hpp:123
@ waterPhaseIdx
Definition: FlowProblem.hpp:136
@ numComponents
Definition: FlowProblem.hpp:116
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: FlowProblem.hpp:101
typename EclThermalLawManager::SolidEnergyLawParams SolidEnergyLawParams
Definition: FlowProblem.hpp:152
bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1453
GetPropType< TypeTag, Properties::BaseProblem > ParentType
Definition: FlowProblem.hpp:98
GetPropType< TypeTag, Properties::EqVector > EqVector
Definition: FlowProblem.hpp:106
unsigned satnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:826
virtual void handleOxygBC(const BCProp::BCFace &, RateVector &) const =0
virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart)=0
virtual void handleUreaBC(const BCProp::BCFace &, RateVector &) const =0
Scalar thermalHalfTransmissibilityOut(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:610
@ enableBrine
Definition: FlowProblem.hpp:119
bool first_step_
Definition: FlowProblem.hpp:1721
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:148
const ThermalConductionLawParams & thermalConductionLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:887
AquiferModel aquiferModel_
Definition: FlowProblem.hpp:1675
GlobalEqVector drift_
Definition: FlowProblem.hpp:1672
bool updateMinPressure_()
Definition: FlowProblem.hpp:1276
std::function< std::vector< double >(const FieldPropsManager &, const std::string &)> fieldPropDoubleOnLeafAssigner_()
Definition: FlowProblem.hpp:1309
Scalar transmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:567
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:145
void updateReferencePorosity_()
Definition: FlowProblem.hpp:1383
Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn, const unsigned globalSpaceIdxOut) const
Definition: FlowProblem.hpp:587
BCData< int > bcindex_
Definition: FlowProblem.hpp:1719
GetPropType< TypeTag, Properties::TracerModel > TracerModel
Definition: FlowProblem.hpp:166
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:129
Scalar rockCompressibility(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:686
bool updateMaxWaterSaturation_()
Definition: FlowProblem.hpp:1245
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: FlowProblem.hpp:163
Scalar maxOilSaturation(unsigned globalDofIdx) const
Returns an element's historic maximum oil phase saturation that was observed during the simulation.
Definition: FlowProblem.hpp:902
bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1291
std::string name() const
The problem name.
Definition: FlowProblem.hpp:857
int episodeIndex() const
Definition: FlowProblem.hpp:280
LhsEval rockCompTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1096
void endIteration()
Called by the simulator after each Newton-Raphson iteration.
Definition: FlowProblem.hpp:388
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:157
FlowProblem(Simulator &simulator)
Definition: FlowProblem.hpp:209
@ enablePolymerMolarWeight
Definition: FlowProblem.hpp:128
@ enableSolvent
Definition: FlowProblem.hpp:130
@ oilPhaseIdx
Definition: FlowProblem.hpp:135
virtual void handleSolventBC(const BCProp::BCFace &, RateVector &) const =0
@ enableDiffusion
Definition: FlowProblem.hpp:120
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: FlowProblem.hpp:105
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: FlowProblem.hpp:146
const Vanguard::TransmissibilityType & eclTransmissibilities() const
Return a reference to the object that handles the "raw" transmissibilities.
Definition: FlowProblem.hpp:633
void source(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:972
virtual void readEquilInitialCondition_()=0
Scalar nextTimeStepSize() const
Propose the size of the next time step to the simulator.
Definition: FlowProblem.hpp:1025
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1630
std::pair< BCType, RateVector > boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
Definition: FlowProblem.hpp:1120
@ enableMICP
Definition: FlowProblem.hpp:126
typename EclMaterialLawManager::MaterialLawParams MaterialLawParams
Definition: FlowProblem.hpp:151
virtual ~FlowProblem()=default
PffGridVector< GridView, Stencil, PffDofData_, DofMapper > pffDofData_
Definition: FlowProblem.hpp:1677
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:778
@ enableFoam
Definition: FlowProblem.hpp:125
const ThermalConductionLawParams & thermalConductionLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:762
Scalar transmissibility(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:509
@ numEq
Definition: FlowProblem.hpp:114
TracerModel tracerModel_
Definition: FlowProblem.hpp:1678
virtual void handlePolymerBC(const BCProp::BCFace &, RateVector &) const =0
@ gasCompIdx
Definition: FlowProblem.hpp:140
const TracerModel & tracerModel() const
Definition: FlowProblem.hpp:637
WellModel wellModel_
Definition: FlowProblem.hpp:1674
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:288
const SolidEnergyLawParams & solidEnergyLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:881
@ dim
Definition: FlowProblem.hpp:110
Scalar getRockCompTransMultVal(std::size_t dofIdx) const
Definition: FlowProblem.hpp:1461
Scalar dofCenterDepth(unsigned globalSpaceIdx) const
Direct indexed acces to the depth of an degree of freedom [m].
Definition: FlowProblem.hpp:677
typename GetProp< TypeTag, Properties::MaterialLaw >::EclMaterialLawManager EclMaterialLawManager
Definition: FlowProblem.hpp:149
LhsEval wellTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Return the well transmissibility multiplier due to rock changes.
Definition: FlowProblem.hpp:1108
const MaterialLawParams & materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
Definition: FlowProblem.hpp:739
std::shared_ptr< const EclMaterialLawManager > materialLawManager() const
Returns the ECL material law manager.
Definition: FlowProblem.hpp:774
unsigned plmixnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:842
const MaterialLawParams & materialLawParams(unsigned globalDofIdx) const
Definition: FlowProblem.hpp:734
Scalar temperature(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:864
Scalar thermalHalfTransmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:623
unsigned miscnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:834
void updateRockCompTransMultVal_()
Definition: FlowProblem.hpp:1612
std::shared_ptr< EclThermalLawManager > thermalLawManager_
Definition: FlowProblem.hpp:1670
Scalar limitNextTimeStepSize_(Scalar dtNext) const
Definition: FlowProblem.hpp:1554
GetPropType< TypeTag, Properties::Stencil > Stencil
Definition: FlowProblem.hpp:103
typename GetProp< TypeTag, Properties::SolidEnergyLaw >::EclThermalLawManager EclThermalLawManager
Definition: FlowProblem.hpp:150
virtual void readExplicitInitialCondition_()=0
virtual void handleMicrBC(const BCProp::BCFace &, RateVector &) const =0
GetPropType< TypeTag, Properties::WellModel > WellModel
Definition: FlowProblem.hpp:159
const DimMatrix & intrinsicPermeability(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:488
bool updateHysteresis_()
Definition: FlowProblem.hpp:1437
void readThermalParameters_()
Definition: FlowProblem.hpp:1367
@ enableEnergy
Definition: FlowProblem.hpp:122
void serializeOp(Serializer &serializer)
Definition: FlowProblem.hpp:1175
@ numPhases
Definition: FlowProblem.hpp:115
Scalar maxPolymerAdsorption(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the max polymer adsorption value.
Definition: FlowProblem.hpp:851
void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:345
Scalar thermalHalfTransmissibilityIn(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:597
void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
Sets an element's maximum oil phase saturation observed during the simulation.
Definition: FlowProblem.hpp:919
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:132
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:749
typename EclThermalLawManager::ThermalConductionLawParams ThermalConductionLawParams
Definition: FlowProblem.hpp:153
Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:547
AquiferModel & mutableAquiferModel()
Definition: FlowProblem.hpp:1013
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:158
Scalar temperature(unsigned globalDofIdx, unsigned) const
Definition: FlowProblem.hpp:873
GetPropType< TypeTag, Properties::AquiferModel > AquiferModel
Definition: FlowProblem.hpp:160
@ waterCompIdx
Definition: FlowProblem.hpp:142
std::shared_ptr< EclMaterialLawManager > materialLawManager_
Definition: FlowProblem.hpp:1669
WellModel & wellModel()
Definition: FlowProblem.hpp:1007
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowProblem.hpp:102
void updateProperty_(const std::string &failureMsg, UpdateFunc func)
Definition: FlowProblem.hpp:1194
bool updateMaxOilSaturation_()
Definition: FlowProblem.hpp:1213
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblem.hpp:180
void updatePffDofData_()
Definition: FlowProblem.hpp:1480
virtual void endEpisode()
Called by the simulator after the end of an episode.
Definition: FlowProblem.hpp:450
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:144
bool nonTrivialBoundaryConditions_
Definition: FlowProblem.hpp:1720
@ enableDispersion
Definition: FlowProblem.hpp:121
@ gasPhaseIdx
Definition: FlowProblem.hpp:134
GetPropType< TypeTag, Properties::Problem > Implementation
Definition: FlowProblem.hpp:99
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1511
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1667
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:963
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:398
Utility::CopyablePtr< DirectionalMobility< TypeTag > > DirectionalMobilityPtr
Definition: FlowProblem.hpp:167
virtual void readInitialCondition_()
Definition: FlowProblem.hpp:1408
virtual void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblem.hpp:930
const MaterialLawParams & materialLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:727
Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:554
@ enableTemperature
Definition: FlowProblem.hpp:131
Scalar diffusivity(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:529
Scalar rockReferencePressure(unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:705
void deserialize(Restarter &res)
This method restores the complete state of the problem and its sub-objects from disk.
Definition: FlowProblem.hpp:254
std::shared_ptr< EclMaterialLawManager > materialLawManager()
Definition: FlowProblem.hpp:810
@ enableConvectiveMixing
Definition: FlowProblem.hpp:118
bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1262
typename GridView::template Codim< 0 >::Entity Element
Definition: FlowProblem.hpp:147
Scalar dofCenterDepth(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the depth of an degree of freedom [m].
Definition: FlowProblem.hpp:665
bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1231
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: FlowProblem.hpp:104
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: FlowProblem.hpp:154
const AquiferModel & aquiferModel() const
Definition: FlowProblem.hpp:1010
int refPressurePhaseIdx_() const
Definition: FlowProblem.hpp:1600
LhsEval rockCompPoroMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the porosity multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1052
Scalar transmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:577
TracerModel & tracerModel()
Definition: FlowProblem.hpp:640
void readMaterialParameters_()
Definition: FlowProblem.hpp:1333
void serialize(Restarter &res)
This method writes the complete state of the problem and its subobjects to disk.
Definition: FlowProblem.hpp:273
MathToolbox< Evaluation > Toolbox
Definition: FlowProblem.hpp:162
Scalar rockReferencePressure(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:696
Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:540
void prefetch(const Element &elem) const
Definition: FlowProblem.hpp:239
A random-access container which stores data attached to a grid's degrees of freedom in a prefetch fri...
Definition: pffgridvector.hh:50
Definition: RelpermDiagnostics.hpp:51
void diagnosis(const EclipseState &eclState, const LevelCartesianIndexMapper &levelCartesianIndexMapper)
This file contains definitions related to directional mobilities.
@ NONE
Definition: DeferredLogger.hpp:46
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: blackoilboundaryratevector.hh:39
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:233
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:224
Definition: FlowProblem.hpp:1682
const std::vector< T > & operator()(FaceDir::DirEnum dir) const
Definition: FlowProblem.hpp:1691
void resize(std::size_t size, T defVal)
Definition: FlowProblem.hpp:1685
std::vector< T > & operator()(FaceDir::DirEnum dir)
Definition: FlowProblem.hpp:1703
std::array< std::vector< T >, 6 > data
Definition: FlowProblem.hpp:1683
Definition: FlowProblem.hpp:1471
ConditionalStorage< enableEnergy, Scalar > thermalHalfTransIn
Definition: FlowProblem.hpp:1472
ConditionalStorage< enableEnergy, Scalar > thermalHalfTransOut
Definition: FlowProblem.hpp:1473
ConditionalStorage< enableDiffusion, Scalar > diffusivity
Definition: FlowProblem.hpp:1474
ConditionalStorage< enableDispersion, Scalar > dispersivity
Definition: FlowProblem.hpp:1475
Scalar transmissibility
Definition: FlowProblem.hpp:1476