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
196 static int handlePositionalParameter(std::function<void(const std::string&,
197 const std::string&)> addKey,
198 std::set<std::string>& seenParams,
199 std::string& errorMsg,
200 int,
201 const char** argv,
202 int paramIdx,
203 int)
204 {
205 return detail::eclPositionalParameter(addKey,
206 seenParams,
207 errorMsg,
208 argv,
209 paramIdx);
210 }
211
215 explicit FlowProblem(Simulator& simulator)
216 : ParentType(simulator)
217 , BaseType(simulator.vanguard().eclState(),
218 simulator.vanguard().schedule(),
219 simulator.vanguard().gridView())
220 , transmissibilities_(simulator.vanguard().eclState(),
221 simulator.vanguard().gridView(),
222 simulator.vanguard().cartesianIndexMapper(),
223 simulator.vanguard().grid(),
224 simulator.vanguard().cellCentroids(),
228 , wellModel_(simulator)
229 , aquiferModel_(simulator)
230 , pffDofData_(simulator.gridView(), this->elementMapper())
231 , tracerModel_(simulator)
232 {
233 if (! Parameters::Get<Parameters::CheckSatfuncConsistency>()) {
234 // User did not enable the "new" saturation function consistency
235 // check module. Run the original checker instead. This is a
236 // temporary measure.
237 RelpermDiagnostics relpermDiagnostics{};
238 relpermDiagnostics.diagnosis(simulator.vanguard().eclState(),
239 simulator.vanguard().levelCartesianIndexMapper());
240 }
241 }
242
243 virtual ~FlowProblem() = default;
244
245 void prefetch(const Element& elem) const
246 { this->pffDofData_.prefetch(elem); }
247
259 template <class Restarter>
260 void deserialize(Restarter& res)
261 {
262 // reload the current episode/report step from the deck
263 this->beginEpisode();
264
265 // deserialize the wells
266 wellModel_.deserialize(res);
267
268 // deserialize the aquifer
269 aquiferModel_.deserialize(res);
270 }
271
278 template <class Restarter>
279 void serialize(Restarter& res)
280 {
281 wellModel_.serialize(res);
282
283 aquiferModel_.serialize(res);
284 }
285
286 int episodeIndex() const
287 {
288 return std::max(this->simulator().episodeIndex(), 0);
289 }
290
294 virtual void beginEpisode()
295 {
296 OPM_TIMEBLOCK(beginEpisode);
297 // Proceed to the next report step
298 auto& simulator = this->simulator();
299 int episodeIdx = simulator.episodeIndex();
300 auto& eclState = simulator.vanguard().eclState();
301 const auto& schedule = simulator.vanguard().schedule();
302 const auto& events = schedule[episodeIdx].events();
303
304 if (episodeIdx >= 0 && events.hasEvent(ScheduleEvents::GEO_MODIFIER)) {
305 // bring the contents of the keywords to the current state of the SCHEDULE
306 // section.
307 //
308 // TODO (?): make grid topology changes possible (depending on what exactly
309 // has changed, the grid may need be re-created which has some serious
310 // implications on e.g., the solution of the simulation.)
311 const auto& miniDeck = schedule[episodeIdx].geo_keywords();
312 const auto& cc = simulator.vanguard().grid().comm();
313 eclState.apply_schedule_keywords( miniDeck );
314 eclBroadcast(cc, eclState.getTransMult() );
315
316 // Re-ordering in case of ALUGrid
317 std::function<unsigned int(unsigned int)> equilGridToGrid = [&simulator](unsigned int i) {
318 return simulator.vanguard().gridEquilIdxToGridIdx(i);
319 };
320
321 // re-compute all quantities which may possibly be affected.
322 using TransUpdateQuantities = typename Vanguard::TransmissibilityType::TransUpdateQuantities;
323 transmissibilities_.update(true, TransUpdateQuantities::All, equilGridToGrid);
324 this->referencePorosity_[1] = this->referencePorosity_[0];
327 this->model().linearizer().updateDiscretizationParameters();
328 }
329
330 bool tuningEvent = this->beginEpisode_(enableExperiments, this->episodeIndex());
331
332 // set up the wells for the next episode.
333 wellModel_.beginEpisode();
334
335 // set up the aquifers for the next episode.
336 aquiferModel_.beginEpisode();
337
338 // set the size of the initial time step of the episode
339 Scalar dt = limitNextTimeStepSize_(simulator.episodeLength());
340 // negative value of initialTimeStepSize_ indicates no active limit from TSINIT or NEXTSTEP
341 if ( (episodeIdx == 0 || tuningEvent) && this->initialTimeStepSize_ > 0)
342 // allow the size of the initial time step to be set via an external parameter
343 // if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
344 dt = std::min(dt, this->initialTimeStepSize_);
345 simulator.setTimeStepSize(dt);
346 }
347
352 {
353 OPM_TIMEBLOCK(beginTimeStep);
354 const int episodeIdx = this->episodeIndex();
355 const int timeStepSize = this->simulator().timeStepSize();
356
358 episodeIdx,
359 this->simulator().timeStepIndex(),
360 this->simulator().startTime(),
361 this->simulator().time(),
362 timeStepSize,
363 this->simulator().endTime());
364
365 // update maximum water saturation and minimum pressure
366 // used when ROCKCOMP is activated
367 // Do not update max RS first step after a restart
368 this->updateExplicitQuantities_(episodeIdx, timeStepSize, first_step_ && (episodeIdx > 0));
369 first_step_ = false;
370
372 this->model().linearizer().updateBoundaryConditionData();
373 }
374
375 wellModel_.beginTimeStep();
376 aquiferModel_.beginTimeStep();
377 tracerModel_.beginTimeStep();
378
379 }
380
385 {
386 OPM_TIMEBLOCK(beginIteration);
387 wellModel_.beginIteration();
388 aquiferModel_.beginIteration();
389 }
390
395 {
396 OPM_TIMEBLOCK(endIteration);
397 wellModel_.endIteration();
398 aquiferModel_.endIteration();
399 }
400
404 virtual void endTimeStep()
405 {
406 OPM_TIMEBLOCK(endTimeStep);
407
408#ifndef NDEBUG
409 if constexpr (getPropValue<TypeTag, Properties::EnableDebuggingChecks>()) {
410 // in debug mode, we don't care about performance, so we check
411 // if the model does the right thing (i.e., the mass change
412 // inside the whole reservoir must be equivalent to the fluxes
413 // over the grid's boundaries plus the source rates specified by
414 // the problem).
415 const int rank = this->simulator().gridView().comm().rank();
416 if (rank == 0) {
417 std::cout << "checking conservativeness of solution\n";
418 }
419
420 this->model().checkConservativeness(/*tolerance=*/-1, /*verbose=*/true);
421 if (rank == 0) {
422 std::cout << "solution is sufficiently conservative\n";
423 }
424 }
425#endif // NDEBUG
426
427 auto& simulator = this->simulator();
428 simulator.setTimeStepIndex(simulator.timeStepIndex()+1);
429
430 this->wellModel_.endTimeStep();
431 this->aquiferModel_.endTimeStep();
432 this->tracerModel_.endTimeStep();
433
434 // Compute flux for output
435 this->model().linearizer().updateFlowsInfo();
436
437 if (this->enableDriftCompensation_) {
438 OPM_TIMEBLOCK(driftCompansation);
439
440 const auto& residual = this->model().linearizer().residual();
441
442 for (unsigned globalDofIdx = 0; globalDofIdx < residual.size(); globalDofIdx ++) {
443 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(globalDofIdx);
444 this->drift_[sfcdofIdx] = residual[sfcdofIdx] * simulator.timeStepSize();
445
446 if constexpr (getPropValue<TypeTag, Properties::UseVolumetricResidual>()) {
447 this->drift_[sfcdofIdx] *= this->model().dofTotalVolume(sfcdofIdx);
448 }
449 }
450 }
451 }
452
456 virtual void endEpisode()
457 {
458 const int episodeIdx = this->episodeIndex();
459
460 this->wellModel_.endEpisode();
461 this->aquiferModel_.endEpisode();
462
463 const auto& schedule = this->simulator().vanguard().schedule();
464
465 // End simulation when completed.
466 if (episodeIdx + 1 >= static_cast<int>(schedule.size()) - 1) {
467 this->simulator().setFinished(true);
468 return;
469 }
470
471 // Otherwise, start next episode (report step).
472 this->simulator().startNextEpisode(schedule.stepLength(episodeIdx + 1));
473 }
474
479 virtual void writeOutput(bool verbose)
480 {
481 OPM_TIMEBLOCK(problemWriteOutput);
482
483 if (Parameters::Get<Parameters::EnableWriteAllSolutions>() ||
484 this->episodeWillBeOver())
485 {
486 // Create VTK output as needed.
487 ParentType::writeOutput(verbose);
488 }
489 }
490
494 template <class Context>
495 const DimMatrix& intrinsicPermeability(const Context& context,
496 unsigned spaceIdx,
497 unsigned timeIdx) const
498 {
499 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
500 return transmissibilities_.permeability(globalSpaceIdx);
501 }
502
509 const DimMatrix& intrinsicPermeability(unsigned globalElemIdx) const
510 { return transmissibilities_.permeability(globalElemIdx); }
511
515 template <class Context>
516 Scalar transmissibility(const Context& context,
517 [[maybe_unused]] unsigned fromDofLocalIdx,
518 unsigned toDofLocalIdx) const
519 {
520 assert(fromDofLocalIdx == 0);
521 return pffDofData_.get(context.element(), toDofLocalIdx).transmissibility;
522 }
523
527 Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
528 {
529 return transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
530 }
531
535 template <class Context>
536 Scalar diffusivity(const Context& context,
537 [[maybe_unused]] unsigned fromDofLocalIdx,
538 unsigned toDofLocalIdx) const
539 {
540 assert(fromDofLocalIdx == 0);
541 return *pffDofData_.get(context.element(), toDofLocalIdx).diffusivity;
542 }
543
547 Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
548 return transmissibilities_.diffusivity(globalCellIn, globalCellOut);
549 }
550
554 Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
555 return transmissibilities_.dispersivity(globalCellIn, globalCellOut);
556 }
557
561 Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx,
562 const unsigned boundaryFaceIdx) const
563 {
564 return transmissibilities_.thermalTransmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
565 }
566
567
568
569
573 template <class Context>
574 Scalar transmissibilityBoundary(const Context& elemCtx,
575 unsigned boundaryFaceIdx) const
576 {
577 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
578 return transmissibilities_.transmissibilityBoundary(elemIdx, boundaryFaceIdx);
579 }
580
584 Scalar transmissibilityBoundary(const unsigned globalSpaceIdx,
585 const unsigned boundaryFaceIdx) const
586 {
587 return transmissibilities_.transmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
588 }
589
590
594 Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn,
595 const unsigned globalSpaceIdxOut) const
596 {
597 return transmissibilities_.thermalHalfTrans(globalSpaceIdxIn,globalSpaceIdxOut);
598 }
599
603 template <class Context>
604 Scalar thermalHalfTransmissibilityIn(const Context& context,
605 unsigned faceIdx,
606 unsigned timeIdx) const
607 {
608 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
609 unsigned toDofLocalIdx = face.exteriorIndex();
610 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransIn;
611 }
612
616 template <class Context>
618 unsigned faceIdx,
619 unsigned timeIdx) const
620 {
621 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
622 unsigned toDofLocalIdx = face.exteriorIndex();
623 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransOut;
624 }
625
629 template <class Context>
631 unsigned boundaryFaceIdx) const
632 {
633 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
634 return transmissibilities_.thermalHalfTransBoundary(elemIdx, boundaryFaceIdx);
635 }
636
640 const typename Vanguard::TransmissibilityType& eclTransmissibilities() const
641 { return transmissibilities_; }
642
643
645 { return tracerModel_; }
646
648 { return tracerModel_; }
649
658 template <class Context>
659 Scalar porosity(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
660 {
661 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
662 return this->porosity(globalSpaceIdx, timeIdx);
663 }
664
671 template <class Context>
672 Scalar dofCenterDepth(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
673 {
674 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
675 return this->dofCenterDepth(globalSpaceIdx);
676 }
677
684 Scalar dofCenterDepth(unsigned globalSpaceIdx) const
685 {
686 return this->simulator().vanguard().cellCenterDepth(globalSpaceIdx);
687 }
688
692 template <class Context>
693 Scalar rockCompressibility(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
694 {
695 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
696 return this->rockCompressibility(globalSpaceIdx);
697 }
698
702 template <class Context>
703 Scalar rockReferencePressure(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
704 {
705 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
706 return rockReferencePressure(globalSpaceIdx);
707 }
708
712 Scalar rockReferencePressure(unsigned globalSpaceIdx) const
713 {
714 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
715 if (rock_config.store()) {
716 return asImp_().initialFluidState(globalSpaceIdx).pressure(refPressurePhaseIdx_());
717 }
718 else {
719 if (this->rockParams_.empty())
720 return 1e5;
721
722 unsigned tableIdx = 0;
723 if (!this->rockTableIdx_.empty()) {
724 tableIdx = this->rockTableIdx_[globalSpaceIdx];
725 }
726 return this->rockParams_[tableIdx].referencePressure;
727 }
728 }
729
733 template <class Context>
734 const MaterialLawParams& materialLawParams(const Context& context,
735 unsigned spaceIdx, unsigned timeIdx) const
736 {
737 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
738 return this->materialLawParams(globalSpaceIdx);
739 }
740
741 const MaterialLawParams& materialLawParams(unsigned globalDofIdx) const
742 {
743 return materialLawManager_->materialLawParams(globalDofIdx);
744 }
745
746 const MaterialLawParams& materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
747 {
748 return materialLawManager_->materialLawParams(globalDofIdx, facedir);
749 }
750
754 template <class Context>
756 solidEnergyLawParams(const Context& context,
757 unsigned spaceIdx,
758 unsigned timeIdx) const
759 {
760 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
761 return thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
762 }
763
767 template <class Context>
769 thermalConductionLawParams(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
770 {
771 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
772 return thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
773 }
774
781 std::shared_ptr<const EclMaterialLawManager> materialLawManager() const
782 { return materialLawManager_; }
783
784 template <class FluidState, class ...Args>
786 std::array<Evaluation,numPhases> &mobility,
788 FluidState &fluidState,
789 unsigned globalSpaceIdx) const
790 {
791 using ContainerT = std::array<Evaluation, numPhases>;
792 OPM_TIMEBLOCK_LOCAL(updateRelperms);
793 {
794 // calculate relative permeabilities. note that we store the result into the
795 // mobility_ class attribute. the division by the phase viscosity happens later.
796 const auto& materialParams = materialLawParams(globalSpaceIdx);
797 MaterialLaw::template relativePermeabilities<ContainerT, FluidState, Args...>(mobility, materialParams, fluidState);
798 Valgrind::CheckDefined(mobility);
799 }
800 if (materialLawManager_->hasDirectionalRelperms()
801 || materialLawManager_->hasDirectionalImbnum())
802 {
803 using Dir = FaceDir::DirEnum;
804 constexpr int ndim = 3;
805 dirMob = std::make_unique<DirectionalMobility<TypeTag>>();
806 Dir facedirs[ndim] = {Dir::XPlus, Dir::YPlus, Dir::ZPlus};
807 for (int i = 0; i<ndim; i++) {
808 const auto& materialParams = materialLawParams(globalSpaceIdx, facedirs[i]);
809 auto& mob_array = dirMob->getArray(i);
810 MaterialLaw::template relativePermeabilities<ContainerT, FluidState, Args...>(mob_array, materialParams, fluidState);
811 }
812 }
813 }
814
818 std::shared_ptr<EclMaterialLawManager> materialLawManager()
819 { return materialLawManager_; }
820
825 template <class Context>
826 unsigned pvtRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
827 { return pvtRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
828
833 template <class Context>
834 unsigned satnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
835 { return this->satnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
836
841 template <class Context>
842 unsigned miscnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
843 { return this->miscnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
844
849 template <class Context>
850 unsigned plmixnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
851 { return this->plmixnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
852
853 // TODO: polymer related might need to go to the blackoil side
858 template <class Context>
859 Scalar maxPolymerAdsorption(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
860 { return this->maxPolymerAdsorption(context.globalSpaceIndex(spaceIdx, timeIdx)); }
861
865 std::string name() const
866 { return this->simulator().vanguard().caseName(); }
867
871 template <class Context>
872 Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
873 {
874 // use the initial temperature of the DOF if temperature is not a primary
875 // variable
876 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
877 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
878 }
879
880
881 Scalar temperature(unsigned globalDofIdx, unsigned /*timeIdx*/) const
882 {
883 // use the initial temperature of the DOF if temperature is not a primary
884 // variable
885 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
886 }
887
889 solidEnergyLawParams(unsigned globalSpaceIdx,
890 unsigned /*timeIdx*/) const
891 {
892 return this->thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
893 }
895 thermalConductionLawParams(unsigned globalSpaceIdx,
896 unsigned /*timeIdx*/)const
897 {
898 return this->thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
899 }
900
910 Scalar maxOilSaturation(unsigned globalDofIdx) const
911 {
912 if (!this->vapparsActive(this->episodeIndex()))
913 return 0.0;
914
915 return this->maxOilSaturation_[globalDofIdx];
916 }
917
927 void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
928 {
929 if (!this->vapparsActive(this->episodeIndex()))
930 return;
931
932 this->maxOilSaturation_[globalDofIdx] = value;
933 }
934
939 {
940 // Calculate all intensive quantities.
941 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/0);
942
943 // We also need the intensive quantities for timeIdx == 1
944 // corresponding to the start of the current timestep, if we
945 // do not use the storage cache, or if we cannot recycle the
946 // first iteration storage.
947 if (!this->model().enableStorageCache() || !this->recycleFirstIterationStorage()) {
948 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/1);
949 }
950
951 // initialize the wells. Note that this needs to be done after initializing the
952 // intrinsic permeabilities and the after applying the initial solution because
953 // the well model uses these...
954 wellModel_.init();
955
956 aquiferModel_.initialSolutionApplied();
957
958 const bool invalidateFromHyst = updateHysteresis_();
959 if (invalidateFromHyst) {
960 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
961 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
962 }
963 }
964
970 template <class Context>
971 void source(RateVector& rate,
972 const Context& context,
973 unsigned spaceIdx,
974 unsigned timeIdx) const
975 {
976 const unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
977 source(rate, globalDofIdx, timeIdx);
978 }
979
980 void source(RateVector& rate,
981 unsigned globalDofIdx,
982 unsigned timeIdx) const
983 {
984 OPM_TIMEBLOCK_LOCAL(eclProblemSource);
985 rate = 0.0;
986
987 // Add well contribution to source here.
988 wellModel_.computeTotalRatesForDof(rate, globalDofIdx);
989
990 // convert the source term from the total mass rate of the
991 // cell to the one per unit of volume as used by the model.
992 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
993 rate[eqIdx] /= this->model().dofTotalVolume(globalDofIdx);
994
995 Valgrind::CheckDefined(rate[eqIdx]);
996 assert(isfinite(rate[eqIdx]));
997 }
998
999 // Add non-well sources.
1000 addToSourceDense(rate, globalDofIdx, timeIdx);
1001 }
1002
1003 virtual void addToSourceDense(RateVector& rate,
1004 unsigned globalDofIdx,
1005 unsigned timeIdx) const = 0;
1006
1012 const WellModel& wellModel() const
1013 { return wellModel_; }
1014
1016 { return wellModel_; }
1017
1019 { return aquiferModel_; }
1020
1022 { return aquiferModel_; }
1023
1026
1034 {
1035 OPM_TIMEBLOCK(nexTimeStepSize);
1036 // allow external code to do the timestepping
1037 if (this->nextTimeStepSize_ > 0.0)
1038 return this->nextTimeStepSize_;
1039
1040 const auto& simulator = this->simulator();
1041 int episodeIdx = simulator.episodeIndex();
1042
1043 // for the initial episode, we use a fixed time step size
1044 if (episodeIdx < 0)
1045 return this->initialTimeStepSize_;
1046
1047 // ask the newton method for a suggestion. This suggestion will be based on how
1048 // well the previous time step converged. After that, apply the runtime time
1049 // stepping constraints.
1050 const auto& newtonMethod = this->model().newtonMethod();
1051 return limitNextTimeStepSize_(newtonMethod.suggestTimeStepSize(simulator.timeStepSize()));
1052 }
1053
1059 template <class LhsEval>
1060 LhsEval rockCompPoroMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1061 {
1062 OPM_TIMEBLOCK_LOCAL(rockCompPoroMultiplier);
1063 if (this->rockCompPoroMult_.empty() && this->rockCompPoroMultWc_.empty())
1064 return 1.0;
1065
1066 unsigned tableIdx = 0;
1067 if (!this->rockTableIdx_.empty())
1068 tableIdx = this->rockTableIdx_[elementIdx];
1069
1070 const auto& fs = intQuants.fluidState();
1071 LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
1072 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1073 if (!this->minRefPressure_.empty())
1074 // The pore space change is irreversible
1075 effectivePressure =
1076 min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
1077 this->minRefPressure_[elementIdx]);
1078
1079 if (!this->overburdenPressure_.empty())
1080 effectivePressure -= this->overburdenPressure_[elementIdx];
1081
1082 if (rock_config.store()) {
1083 effectivePressure -= asImp_().initialFluidState(elementIdx).pressure(refPressurePhaseIdx_());
1084 }
1085
1086 if (!this->rockCompPoroMult_.empty()) {
1087 return this->rockCompPoroMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1088 }
1089
1090 // water compaction
1091 assert(!this->rockCompPoroMultWc_.empty());
1092 LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1093 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1094
1095 return this->rockCompPoroMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1096 }
1097
1103 template <class LhsEval>
1104 LhsEval rockCompTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1105 {
1106 const bool implicit = !this->explicitRockCompaction_;
1107 return implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx)
1108 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1109 }
1110
1111
1115 template <class LhsEval>
1116 LhsEval wellTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1117 {
1118 OPM_TIMEBLOCK_LOCAL(wellTransMultiplier);
1119
1120 const bool implicit = !this->explicitRockCompaction_;
1121 double trans_mult = implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<double>(intQuants, elementIdx)
1122 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1123 trans_mult *= this->simulator().problem().template permFactTransMultiplier<double>(intQuants, elementIdx);
1124
1125 return trans_mult;
1126 }
1127
1128 std::pair<BCType, RateVector> boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
1129 {
1130 OPM_TIMEBLOCK_LOCAL(boundaryCondition);
1132 return { BCType::NONE, RateVector(0.0) };
1133 }
1134 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
1135 const auto& schedule = this->simulator().vanguard().schedule();
1136 if (bcindex_(dir)[globalSpaceIdx] == 0) {
1137 return { BCType::NONE, RateVector(0.0) };
1138 }
1139 if (schedule[this->episodeIndex()].bcprop.size() == 0) {
1140 return { BCType::NONE, RateVector(0.0) };
1141 }
1142 const auto& bc = schedule[this->episodeIndex()].bcprop[bcindex_(dir)[globalSpaceIdx]];
1143 if (bc.bctype!=BCType::RATE) {
1144 return { bc.bctype, RateVector(0.0) };
1145 }
1146
1147 RateVector rate = 0.0;
1148 switch (bc.component) {
1149 case BCComponent::OIL:
1150 rate[FluidSystem::canonicalToActiveCompIdx(oilCompIdx)] = bc.rate;
1151 break;
1152 case BCComponent::GAS:
1153 rate[FluidSystem::canonicalToActiveCompIdx(gasCompIdx)] = bc.rate;
1154 break;
1155 case BCComponent::WATER:
1156 rate[FluidSystem::canonicalToActiveCompIdx(waterCompIdx)] = bc.rate;
1157 break;
1158 case BCComponent::SOLVENT:
1159 this->handleSolventBC(bc, rate);
1160 break;
1161 case BCComponent::POLYMER:
1162 this->handlePolymerBC(bc, rate);
1163 break;
1164 case BCComponent::MICR:
1165 this->handleMicrBC(bc, rate);
1166 break;
1167 case BCComponent::OXYG:
1168 this->handleOxygBC(bc, rate);
1169 break;
1170 case BCComponent::UREA:
1171 this->handleUreaBC(bc, rate);
1172 break;
1173 case BCComponent::NONE:
1174 throw std::logic_error("you need to specify the component when RATE type is set in BC");
1175 break;
1176 }
1177 //TODO add support for enthalpy rate
1178 return {bc.bctype, rate};
1179 }
1180
1181
1182 template<class Serializer>
1183 void serializeOp(Serializer& serializer)
1184 {
1185 serializer(static_cast<BaseType&>(*this));
1186 serializer(drift_);
1187 serializer(wellModel_);
1188 serializer(aquiferModel_);
1189 serializer(tracerModel_);
1190 serializer(*materialLawManager_);
1191 }
1192
1193private:
1194 Implementation& asImp_()
1195 { return *static_cast<Implementation *>(this); }
1196
1197 const Implementation& asImp_() const
1198 { return *static_cast<const Implementation *>(this); }
1199
1200protected:
1201 template<class UpdateFunc>
1202 void updateProperty_(const std::string& failureMsg,
1203 UpdateFunc func)
1204 {
1205 OPM_TIMEBLOCK(updateProperty);
1206 const auto& model = this->simulator().model();
1207 const auto& primaryVars = model.solution(/*timeIdx*/0);
1208 const auto& vanguard = this->simulator().vanguard();
1209 std::size_t numGridDof = primaryVars.size();
1211#ifdef _OPENMP
1212#pragma omp parallel for
1213#endif
1214 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
1215 const auto& iq = *model.cachedIntensiveQuantities(dofIdx, /*timeIdx=*/ 0);
1216 func(dofIdx, iq);
1217 }
1218 OPM_END_PARALLEL_TRY_CATCH(failureMsg, vanguard.grid().comm());
1219 }
1220
1222 {
1223 OPM_TIMEBLOCK(updateMaxOilSaturation);
1224 int episodeIdx = this->episodeIndex();
1225
1226 // we use VAPPARS
1227 if (this->vapparsActive(episodeIdx)) {
1228 this->updateProperty_("FlowProblem::updateMaxOilSaturation_() failed:",
1229 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1230 {
1231 this->updateMaxOilSaturation_(compressedDofIdx,iq);
1232 });
1233 return true;
1234 }
1235
1236 return false;
1237 }
1238
1239 bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1240 {
1241 OPM_TIMEBLOCK_LOCAL(updateMaxOilSaturation);
1242 const auto& fs = iq.fluidState();
1243 const Scalar So = decay<Scalar>(fs.saturation(refPressurePhaseIdx_()));
1244 auto& mos = this->maxOilSaturation_;
1245 if(mos[compressedDofIdx] < So){
1246 mos[compressedDofIdx] = So;
1247 return true;
1248 }else{
1249 return false;
1250 }
1251 }
1252
1254 {
1255 OPM_TIMEBLOCK(updateMaxWaterSaturation);
1256 // water compaction is activated in ROCKCOMP
1257 if (this->maxWaterSaturation_.empty())
1258 return false;
1259
1260 this->maxWaterSaturation_[/*timeIdx=*/1] = this->maxWaterSaturation_[/*timeIdx=*/0];
1261 this->updateProperty_("FlowProblem::updateMaxWaterSaturation_() failed:",
1262 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1263 {
1264 this->updateMaxWaterSaturation_(compressedDofIdx,iq);
1265 });
1266 return true;
1267 }
1268
1269
1270 bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1271 {
1272 OPM_TIMEBLOCK_LOCAL(updateMaxWaterSaturation);
1273 const auto& fs = iq.fluidState();
1274 const Scalar Sw = decay<Scalar>(fs.saturation(waterPhaseIdx));
1275 auto& mow = this->maxWaterSaturation_;
1276 if(mow[compressedDofIdx]< Sw){
1277 mow[compressedDofIdx] = Sw;
1278 return true;
1279 }else{
1280 return false;
1281 }
1282 }
1283
1285 {
1286 OPM_TIMEBLOCK(updateMinPressure);
1287 // IRREVERS option is used in ROCKCOMP
1288 if (this->minRefPressure_.empty())
1289 return false;
1290
1291 this->updateProperty_("FlowProblem::updateMinPressure_() failed:",
1292 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1293 {
1294 this->updateMinPressure_(compressedDofIdx,iq);
1295 });
1296 return true;
1297 }
1298
1299 bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities& iq){
1300 OPM_TIMEBLOCK_LOCAL(updateMinPressure);
1301 const auto& fs = iq.fluidState();
1302 const Scalar min_pressure = getValue(fs.pressure(refPressurePhaseIdx_()));
1303 auto& min_pressures = this->minRefPressure_;
1304 if(min_pressures[compressedDofIdx]> min_pressure){
1305 min_pressures[compressedDofIdx] = min_pressure;
1306 return true;
1307 }else{
1308 return false;
1309 }
1310 }
1311
1312 // \brief Function to assign field properties of type double, on the leaf grid view.
1313 //
1314 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1315 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1316 std::function<std::vector<double>(const FieldPropsManager&, const std::string&)>
1318 {
1319 const auto& lookup = this->lookUpData_;
1320 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString)
1321 {
1322 return lookup.assignFieldPropsDoubleOnLeaf(fieldPropManager, propString);
1323 };
1324 }
1325
1326 // \brief Function to assign field properties of type int, unsigned int, ..., on the leaf grid view.
1327 //
1328 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1329 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1330 template<typename IntType>
1331 std::function<std::vector<IntType>(const FieldPropsManager&, const std::string&, bool)>
1333 {
1334 const auto& lookup = this->lookUpData_;
1335 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString, bool needsTranslation)
1336 {
1337 return lookup.template assignFieldPropsIntOnLeaf<IntType>(fieldPropManager, propString, needsTranslation);
1338 };
1339 }
1340
1342 {
1343 OPM_TIMEBLOCK(readMaterialParameters);
1344 const auto& simulator = this->simulator();
1345 const auto& vanguard = simulator.vanguard();
1346 const auto& eclState = vanguard.eclState();
1347
1348 // the PVT and saturation region numbers
1350 this->updatePvtnum_();
1351 this->updateSatnum_();
1352
1353 // the MISC region numbers (solvent model)
1354 this->updateMiscnum_();
1355 // the PLMIX region numbers (polymer model)
1356 this->updatePlmixnum_();
1357
1358 OPM_END_PARALLEL_TRY_CATCH("Invalid region numbers: ", vanguard.gridView().comm());
1360 // porosity
1362 this->referencePorosity_[1] = this->referencePorosity_[0];
1364
1366 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1367 materialLawManager_ = std::make_shared<EclMaterialLawManager>();
1368 materialLawManager_->initFromState(eclState);
1369 materialLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1370 this-> template fieldPropIntTypeOnLeafAssigner_<int>(),
1373 }
1374
1376 {
1377 if constexpr (enableEnergy)
1378 {
1379 const auto& simulator = this->simulator();
1380 const auto& vanguard = simulator.vanguard();
1381 const auto& eclState = vanguard.eclState();
1382
1383 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1384 thermalLawManager_ = std::make_shared<EclThermalLawManager>();
1385 thermalLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1387 this-> template fieldPropIntTypeOnLeafAssigner_<unsigned int>());
1388 }
1389 }
1390
1392 {
1393 const auto& simulator = this->simulator();
1394 const auto& vanguard = simulator.vanguard();
1395 const auto& eclState = vanguard.eclState();
1396
1397 std::size_t numDof = this->model().numGridDof();
1398
1399 this->referencePorosity_[/*timeIdx=*/0].resize(numDof);
1400
1401 const auto& fp = eclState.fieldProps();
1402 const std::vector<double> porvData = this -> fieldPropDoubleOnLeafAssigner_()(fp, "PORV");
1403 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1404 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(dofIdx);
1405 Scalar poreVolume = porvData[dofIdx];
1406
1407 // we define the porosity as the accumulated pore volume divided by the
1408 // geometric volume of the element. Note that -- in pathetic cases -- it can
1409 // be larger than 1.0!
1410 Scalar dofVolume = simulator.model().dofTotalVolume(sfcdofIdx);
1411 assert(dofVolume > 0.0);
1412 this->referencePorosity_[/*timeIdx=*/0][sfcdofIdx] = poreVolume/dofVolume;
1413 }
1414 }
1415
1417 {
1418 // TODO: whether we should move this to FlowProblemBlackoil
1419 const auto& simulator = this->simulator();
1420 const auto& vanguard = simulator.vanguard();
1421 const auto& eclState = vanguard.eclState();
1422
1423 if (eclState.getInitConfig().hasEquil())
1425 else
1427
1428 //initialize min/max values
1429 std::size_t numElems = this->model().numGridDof();
1430 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1431 const auto& fs = asImp_().initialFluidStates()[elemIdx];
1432 if (!this->maxWaterSaturation_.empty() && waterPhaseIdx > -1)
1433 this->maxWaterSaturation_[elemIdx] = std::max(this->maxWaterSaturation_[elemIdx], fs.saturation(waterPhaseIdx));
1434 if (!this->maxOilSaturation_.empty() && oilPhaseIdx > -1)
1435 this->maxOilSaturation_[elemIdx] = std::max(this->maxOilSaturation_[elemIdx], fs.saturation(oilPhaseIdx));
1436 if (!this->minRefPressure_.empty() && refPressurePhaseIdx_() > -1)
1437 this->minRefPressure_[elemIdx] = std::min(this->minRefPressure_[elemIdx], fs.pressure(refPressurePhaseIdx_()));
1438 }
1439 }
1440
1441 virtual void readEquilInitialCondition_() = 0;
1443
1444 // update the hysteresis parameters of the material laws for the whole grid
1446 {
1447 if (!materialLawManager_->enableHysteresis())
1448 return false;
1449
1450 // we need to update the hysteresis data for _all_ elements (i.e., not just the
1451 // interior ones) to avoid desynchronization of the processes in the parallel case!
1452 this->updateProperty_("FlowProblem::updateHysteresis_() failed:",
1453 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1454 {
1455 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1456 });
1457 return true;
1458 }
1459
1460
1461 bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1462 {
1463 OPM_TIMEBLOCK_LOCAL(updateHysteresis_);
1464 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1465 //TODO change materials to give a bool
1466 return true;
1467 }
1468
1469 Scalar getRockCompTransMultVal(std::size_t dofIdx) const
1470 {
1471 if (this->rockCompTransMultVal_.empty())
1472 return 1.0;
1473
1474 return this->rockCompTransMultVal_[dofIdx];
1475 }
1476
1477protected:
1479 {
1480 ConditionalStorage<enableEnergy, Scalar> thermalHalfTransIn;
1481 ConditionalStorage<enableEnergy, Scalar> thermalHalfTransOut;
1482 ConditionalStorage<enableDiffusion, Scalar> diffusivity;
1483 ConditionalStorage<enableDispersion, Scalar> dispersivity;
1485 };
1486
1487 // update the prefetch friendly data object
1489 {
1490 const auto& distFn =
1491 [this](PffDofData_& dofData,
1492 const Stencil& stencil,
1493 unsigned localDofIdx)
1494 -> void
1495 {
1496 const auto& elementMapper = this->model().elementMapper();
1497
1498 unsigned globalElemIdx = elementMapper.index(stencil.entity(localDofIdx));
1499 if (localDofIdx != 0) {
1500 unsigned globalCenterElemIdx = elementMapper.index(stencil.entity(/*dofIdx=*/0));
1501 dofData.transmissibility = transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
1502
1503 if constexpr (enableEnergy) {
1504 *dofData.thermalHalfTransIn = transmissibilities_.thermalHalfTrans(globalCenterElemIdx, globalElemIdx);
1505 *dofData.thermalHalfTransOut = transmissibilities_.thermalHalfTrans(globalElemIdx, globalCenterElemIdx);
1506 }
1507 if constexpr (enableDiffusion)
1508 *dofData.diffusivity = transmissibilities_.diffusivity(globalCenterElemIdx, globalElemIdx);
1509 if (enableDispersion)
1510 dofData.dispersivity = transmissibilities_.dispersivity(globalCenterElemIdx, globalElemIdx);
1511 }
1512 };
1513
1514 pffDofData_.update(distFn);
1515 }
1516
1517 virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart) = 0;
1518
1520 {
1521 const auto& simulator = this->simulator();
1522 const auto& vanguard = simulator.vanguard();
1523 const auto& bcconfig = vanguard.eclState().getSimulationConfig().bcconfig();
1524 if (bcconfig.size() > 0) {
1526
1527 std::size_t numCartDof = vanguard.cartesianSize();
1528 unsigned numElems = vanguard.gridView().size(/*codim=*/0);
1529 std::vector<int> cartesianToCompressedElemIdx(numCartDof, -1);
1530
1531 for (unsigned elemIdx = 0; elemIdx < numElems; ++elemIdx)
1532 cartesianToCompressedElemIdx[vanguard.cartesianIndex(elemIdx)] = elemIdx;
1533
1534 bcindex_.resize(numElems, 0);
1535 auto loopAndApply = [&cartesianToCompressedElemIdx,
1536 &vanguard](const auto& bcface,
1537 auto apply)
1538 {
1539 for (int i = bcface.i1; i <= bcface.i2; ++i) {
1540 for (int j = bcface.j1; j <= bcface.j2; ++j) {
1541 for (int k = bcface.k1; k <= bcface.k2; ++k) {
1542 std::array<int, 3> tmp = {i,j,k};
1543 auto elemIdx = cartesianToCompressedElemIdx[vanguard.cartesianIndex(tmp)];
1544 if (elemIdx >= 0)
1545 apply(elemIdx);
1546 }
1547 }
1548 }
1549 };
1550 for (const auto& bcface : bcconfig) {
1551 std::vector<int>& data = bcindex_(bcface.dir);
1552 const int index = bcface.index;
1553 loopAndApply(bcface,
1554 [&data,index](int elemIdx)
1555 { data[elemIdx] = index; });
1556 }
1557 }
1558 }
1559
1560 // this method applies the runtime constraints specified via the deck and/or command
1561 // line parameters for the size of the next time step.
1563 {
1564 if constexpr (enableExperiments) {
1565 const auto& simulator = this->simulator();
1566 const auto& schedule = simulator.vanguard().schedule();
1567 int episodeIdx = simulator.episodeIndex();
1568
1569 // first thing in the morning, limit the time step size to the maximum size
1570 Scalar maxTimeStepSize = Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60;
1571 int reportStepIdx = std::max(episodeIdx, 0);
1572 if (this->enableTuning_) {
1573 const auto& tuning = schedule[reportStepIdx].tuning();
1574 maxTimeStepSize = tuning.TSMAXZ;
1575 }
1576
1577 dtNext = std::min(dtNext, maxTimeStepSize);
1578
1579 Scalar remainingEpisodeTime =
1580 simulator.episodeStartTime() + simulator.episodeLength()
1581 - (simulator.startTime() + simulator.time());
1582 assert(remainingEpisodeTime >= 0.0);
1583
1584 // if we would have a small amount of time left over in the current episode, make
1585 // two equal time steps instead of a big and a small one
1586 if (remainingEpisodeTime/2.0 < dtNext && dtNext < remainingEpisodeTime*(1.0 - 1e-5))
1587 // note: limiting to the maximum time step size here is probably not strictly
1588 // necessary, but it should not hurt and is more fool-proof
1589 dtNext = std::min(maxTimeStepSize, remainingEpisodeTime/2.0);
1590
1591 if (simulator.episodeStarts()) {
1592 // if a well event occurred, respect the limit for the maximum time step after
1593 // that, too
1594 const auto& events = simulator.vanguard().schedule()[reportStepIdx].events();
1595 bool wellEventOccured =
1596 events.hasEvent(ScheduleEvents::NEW_WELL)
1597 || events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE)
1598 || events.hasEvent(ScheduleEvents::INJECTION_UPDATE)
1599 || events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
1600 if (episodeIdx >= 0 && wellEventOccured && this->maxTimeStepAfterWellEvent_ > 0)
1601 dtNext = std::min(dtNext, this->maxTimeStepAfterWellEvent_);
1602 }
1603 }
1604
1605 return dtNext;
1606 }
1607
1609 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1610 return oilPhaseIdx;
1611 }
1612 else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
1613 return gasPhaseIdx;
1614 }
1615 else {
1616 return waterPhaseIdx;
1617 }
1618 }
1619
1621 {
1622 const auto& model = this->simulator().model();
1623 std::size_t numGridDof = this->model().numGridDof();
1624 this->rockCompTransMultVal_.resize(numGridDof, 1.0);
1625 for (std::size_t elementIdx = 0; elementIdx < numGridDof; ++elementIdx) {
1626 const auto& iq = *model.cachedIntensiveQuantities(elementIdx, /*timeIdx=*/ 0);
1627 Scalar trans_mult = computeRockCompTransMultiplier_<Scalar>(iq, elementIdx);
1628 this->rockCompTransMultVal_[elementIdx] = trans_mult;
1629 }
1630 }
1631
1637 template <class LhsEval>
1638 LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1639 {
1640 OPM_TIMEBLOCK_LOCAL(computeRockCompTransMultiplier);
1641 if (this->rockCompTransMult_.empty() && this->rockCompTransMultWc_.empty())
1642 return 1.0;
1643
1644 unsigned tableIdx = 0;
1645 if (!this->rockTableIdx_.empty())
1646 tableIdx = this->rockTableIdx_[elementIdx];
1647
1648 const auto& fs = intQuants.fluidState();
1649 LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
1650 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1651 if (!this->minRefPressure_.empty())
1652 // The pore space change is irreversible
1653 effectivePressure =
1654 min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
1655 this->minRefPressure_[elementIdx]);
1656
1657 if (!this->overburdenPressure_.empty())
1658 effectivePressure -= this->overburdenPressure_[elementIdx];
1659
1660 if (rock_config.store()) {
1661 effectivePressure -= asImp_().initialFluidState(elementIdx).pressure(refPressurePhaseIdx_());
1662 }
1663
1664 if (!this->rockCompTransMult_.empty())
1665 return this->rockCompTransMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1666
1667 // water compaction
1668 assert(!this->rockCompTransMultWc_.empty());
1669 LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1670 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1671
1672 return this->rockCompTransMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1673 }
1674
1675 typename Vanguard::TransmissibilityType transmissibilities_;
1676
1677 std::shared_ptr<EclMaterialLawManager> materialLawManager_;
1678 std::shared_ptr<EclThermalLawManager> thermalLawManager_;
1679
1681
1684
1687
1688 template<class T>
1689 struct BCData
1690 {
1691 std::array<std::vector<T>,6> data;
1692
1693 void resize(std::size_t size, T defVal)
1694 {
1695 for (auto& d : data)
1696 d.resize(size, defVal);
1697 }
1698
1699 const std::vector<T>& operator()(FaceDir::DirEnum dir) const
1700 {
1701 if (dir == FaceDir::DirEnum::Unknown)
1702 throw std::runtime_error("Tried to access BC data for the 'Unknown' direction");
1703 int idx = 0;
1704 int div = static_cast<int>(dir);
1705 while ((div /= 2) >= 1)
1706 ++idx;
1707 assert(idx >= 0 && idx <= 5);
1708 return data[idx];
1709 }
1710
1711 std::vector<T>& operator()(FaceDir::DirEnum dir)
1712 {
1713 return const_cast<std::vector<T>&>(std::as_const(*this)(dir));
1714 }
1715 };
1716
1717 virtual void handleSolventBC(const BCProp::BCFace&, RateVector&) const = 0;
1718
1719 virtual void handlePolymerBC(const BCProp::BCFace&, RateVector&) const = 0;
1720
1721 virtual void handleMicrBC(const BCProp::BCFace&, RateVector&) const = 0;
1722
1723 virtual void handleOxygBC(const BCProp::BCFace&, RateVector&) const = 0;
1724
1725 virtual void handleUreaBC(const BCProp::BCFace&, RateVector&) const = 0;
1726
1729 bool first_step_ = true;
1730
1733 virtual bool episodeWillBeOver() const
1734 {
1735 return this->simulator().episodeWillBeOver();
1736 }
1737};
1738
1739} // namespace Opm
1740
1741#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:736
unsigned pvtRegionIndex(unsigned elemIdx) const
Returns the index the relevant PVT region given a cell index.
Definition: FlowGenericProblem_impl.hpp:695
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:332
Scalar rockCompressibility(unsigned globalSpaceIdx) const
Definition: FlowGenericProblem_impl.hpp:317
unsigned miscnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant MISC region given a cell index.
Definition: FlowGenericProblem_impl.hpp:715
unsigned satnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant saturation function region given a cell index.
Definition: FlowGenericProblem_impl.hpp:705
void beginTimeStep_(bool enableExperiments, int episodeIdx, int timeStepIndex, Scalar startTime, Scalar time, Scalar timeStepSize, Scalar endTime)
Definition: FlowGenericProblem_impl.hpp:451
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:725
std::array< std::vector< Scalar >, 2 > referencePorosity_
Definition: FlowGenericProblem.hpp:319
bool beginEpisode_(bool enableExperiments, int episodeIdx)
Definition: FlowGenericProblem_impl.hpp:414
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:114
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
virtual bool episodeWillBeOver() const
Definition: FlowProblem.hpp:1733
const WellModel & wellModel() const
Returns a reference to the ECL well manager used by the problem.
Definition: FlowProblem.hpp:1012
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:527
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:196
@ dim
Definition: FlowProblem.hpp:110
bool nonTrivialBoundaryConditions() const
Definition: FlowProblem.hpp:1024
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:479
std::function< std::vector< IntType >(const FieldPropsManager &, const std::string &, bool)> fieldPropIntTypeOnLeafAssigner_()
Definition: FlowProblem.hpp:1332
const DimMatrix & intrinsicPermeability(unsigned globalElemIdx) const
This method returns the intrinsic permeability tensor given a global element index.
Definition: FlowProblem.hpp:509
@ gasCompIdx
Definition: FlowProblem.hpp:140
@ enableMICP
Definition: FlowProblem.hpp:126
@ enableExperiments
Definition: FlowProblem.hpp:123
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:826
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:659
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition: FlowProblem.hpp:384
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:107
GetPropType< TypeTag, Properties::DofMapper > DofMapper
Definition: FlowProblem.hpp:155
@ enableEnergy
Definition: FlowProblem.hpp:122
@ oilPhaseIdx
Definition: FlowProblem.hpp:135
@ numPhases
Definition: FlowProblem.hpp:115
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:1461
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:834
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:617
@ enableConvectiveMixing
Definition: FlowProblem.hpp:118
bool first_step_
Definition: FlowProblem.hpp:1729
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:148
const ThermalConductionLawParams & thermalConductionLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:895
AquiferModel aquiferModel_
Definition: FlowProblem.hpp:1683
GlobalEqVector drift_
Definition: FlowProblem.hpp:1680
bool updateMinPressure_()
Definition: FlowProblem.hpp:1284
std::function< std::vector< double >(const FieldPropsManager &, const std::string &)> fieldPropDoubleOnLeafAssigner_()
Definition: FlowProblem.hpp:1317
Scalar transmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:574
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:145
void updateReferencePorosity_()
Definition: FlowProblem.hpp:1391
Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn, const unsigned globalSpaceIdxOut) const
Definition: FlowProblem.hpp:594
BCData< int > bcindex_
Definition: FlowProblem.hpp:1727
GetPropType< TypeTag, Properties::TracerModel > TracerModel
Definition: FlowProblem.hpp:166
@ enablePolymerMolarWeight
Definition: FlowProblem.hpp:128
Scalar rockCompressibility(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:693
bool updateMaxWaterSaturation_()
Definition: FlowProblem.hpp:1253
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:910
bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1299
std::string name() const
The problem name.
Definition: FlowProblem.hpp:865
int episodeIndex() const
Definition: FlowProblem.hpp:286
LhsEval rockCompTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1104
void endIteration()
Called by the simulator after each Newton-Raphson iteration.
Definition: FlowProblem.hpp:394
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:157
FlowProblem(Simulator &simulator)
Definition: FlowProblem.hpp:215
@ enablePolymer
Definition: FlowProblem.hpp:127
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:129
@ gasPhaseIdx
Definition: FlowProblem.hpp:134
virtual void handleSolventBC(const BCProp::BCFace &, RateVector &) const =0
@ enableBrine
Definition: FlowProblem.hpp:119
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:640
void source(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:980
virtual void readEquilInitialCondition_()=0
Scalar nextTimeStepSize() const
Propose the size of the next time step to the simulator.
Definition: FlowProblem.hpp:1033
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1638
std::pair< BCType, RateVector > boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
Definition: FlowProblem.hpp:1128
@ enableFoam
Definition: FlowProblem.hpp:125
typename EclMaterialLawManager::MaterialLawParams MaterialLawParams
Definition: FlowProblem.hpp:151
virtual ~FlowProblem()=default
PffGridVector< GridView, Stencil, PffDofData_, DofMapper > pffDofData_
Definition: FlowProblem.hpp:1685
@ enableExtbo
Definition: FlowProblem.hpp:124
const ThermalConductionLawParams & thermalConductionLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:769
Scalar transmissibility(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:516
@ dimWorld
Definition: FlowProblem.hpp:111
TracerModel tracerModel_
Definition: FlowProblem.hpp:1686
virtual void handlePolymerBC(const BCProp::BCFace &, RateVector &) const =0
@ waterPhaseIdx
Definition: FlowProblem.hpp:136
const TracerModel & tracerModel() const
Definition: FlowProblem.hpp:644
WellModel wellModel_
Definition: FlowProblem.hpp:1682
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:294
const SolidEnergyLawParams & solidEnergyLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:889
Scalar getRockCompTransMultVal(std::size_t dofIdx) const
Definition: FlowProblem.hpp:1469
Scalar dofCenterDepth(unsigned globalSpaceIdx) const
Direct indexed acces to the depth of an degree of freedom [m].
Definition: FlowProblem.hpp:684
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:1116
const MaterialLawParams & materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
Definition: FlowProblem.hpp:746
std::shared_ptr< const EclMaterialLawManager > materialLawManager() const
Returns the ECL material law manager.
Definition: FlowProblem.hpp:781
unsigned plmixnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:850
const MaterialLawParams & materialLawParams(unsigned globalDofIdx) const
Definition: FlowProblem.hpp:741
Scalar temperature(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:872
Scalar thermalHalfTransmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:630
unsigned miscnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:842
void updateRockCompTransMultVal_()
Definition: FlowProblem.hpp:1620
std::shared_ptr< EclThermalLawManager > thermalLawManager_
Definition: FlowProblem.hpp:1678
Scalar limitNextTimeStepSize_(Scalar dtNext) const
Definition: FlowProblem.hpp:1562
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:495
bool updateHysteresis_()
Definition: FlowProblem.hpp:1445
void readThermalParameters_()
Definition: FlowProblem.hpp:1375
@ enableDispersion
Definition: FlowProblem.hpp:121
void serializeOp(Serializer &serializer)
Definition: FlowProblem.hpp:1183
@ numEq
Definition: FlowProblem.hpp:114
Scalar maxPolymerAdsorption(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the max polymer adsorption value.
Definition: FlowProblem.hpp:859
void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:351
Scalar thermalHalfTransmissibilityIn(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:604
void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
Sets an element's maximum oil phase saturation observed during the simulation.
Definition: FlowProblem.hpp:927
@ enableTemperature
Definition: FlowProblem.hpp:131
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:756
typename EclThermalLawManager::ThermalConductionLawParams ThermalConductionLawParams
Definition: FlowProblem.hpp:153
Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:554
AquiferModel & mutableAquiferModel()
Definition: FlowProblem.hpp:1021
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:158
Scalar temperature(unsigned globalDofIdx, unsigned) const
Definition: FlowProblem.hpp:881
GetPropType< TypeTag, Properties::AquiferModel > AquiferModel
Definition: FlowProblem.hpp:160
@ oilCompIdx
Definition: FlowProblem.hpp:141
std::shared_ptr< EclMaterialLawManager > materialLawManager_
Definition: FlowProblem.hpp:1677
WellModel & wellModel()
Definition: FlowProblem.hpp:1015
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowProblem.hpp:102
void updateProperty_(const std::string &failureMsg, UpdateFunc func)
Definition: FlowProblem.hpp:1202
bool updateMaxOilSaturation_()
Definition: FlowProblem.hpp:1221
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblem.hpp:180
void updatePffDofData_()
Definition: FlowProblem.hpp:1488
virtual void endEpisode()
Called by the simulator after the end of an episode.
Definition: FlowProblem.hpp:456
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:144
bool nonTrivialBoundaryConditions_
Definition: FlowProblem.hpp:1728
@ enableDiffusion
Definition: FlowProblem.hpp:120
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:132
GetPropType< TypeTag, Properties::Problem > Implementation
Definition: FlowProblem.hpp:99
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1519
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:785
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1675
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:971
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:404
Utility::CopyablePtr< DirectionalMobility< TypeTag > > DirectionalMobilityPtr
Definition: FlowProblem.hpp:167
virtual void readInitialCondition_()
Definition: FlowProblem.hpp:1416
virtual void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblem.hpp:938
const MaterialLawParams & materialLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:734
Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:561
@ enableSolvent
Definition: FlowProblem.hpp:130
Scalar diffusivity(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:536
Scalar rockReferencePressure(unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:712
@ waterCompIdx
Definition: FlowProblem.hpp:142
void deserialize(Restarter &res)
This method restores the complete state of the problem and its sub-objects from disk.
Definition: FlowProblem.hpp:260
std::shared_ptr< EclMaterialLawManager > materialLawManager()
Definition: FlowProblem.hpp:818
@ numComponents
Definition: FlowProblem.hpp:116
bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1270
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:672
bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1239
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:1018
int refPressurePhaseIdx_() const
Definition: FlowProblem.hpp:1608
LhsEval rockCompPoroMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the porosity multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1060
Scalar transmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:584
TracerModel & tracerModel()
Definition: FlowProblem.hpp:647
void readMaterialParameters_()
Definition: FlowProblem.hpp:1341
void serialize(Restarter &res)
This method writes the complete state of the problem and its subobjects to disk.
Definition: FlowProblem.hpp:279
MathToolbox< Evaluation > Toolbox
Definition: FlowProblem.hpp:162
Scalar rockReferencePressure(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:703
Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:547
void prefetch(const Element &elem) const
Definition: FlowProblem.hpp:245
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:1690
const std::vector< T > & operator()(FaceDir::DirEnum dir) const
Definition: FlowProblem.hpp:1699
void resize(std::size_t size, T defVal)
Definition: FlowProblem.hpp:1693
std::vector< T > & operator()(FaceDir::DirEnum dir)
Definition: FlowProblem.hpp:1711
std::array< std::vector< T >, 6 > data
Definition: FlowProblem.hpp:1691
Definition: FlowProblem.hpp:1479
ConditionalStorage< enableEnergy, Scalar > thermalHalfTransIn
Definition: FlowProblem.hpp:1480
ConditionalStorage< enableEnergy, Scalar > thermalHalfTransOut
Definition: FlowProblem.hpp:1481
ConditionalStorage< enableDiffusion, Scalar > diffusivity
Definition: FlowProblem.hpp:1482
ConditionalStorage< enableDispersion, Scalar > dispersivity
Definition: FlowProblem.hpp:1483
Scalar transmissibility
Definition: FlowProblem.hpp:1484