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
68
72
73#include <opm/utility/CopyablePtr.hpp>
74
75#include <algorithm>
76#include <cstddef>
77#include <functional>
78#include <set>
79#include <stdexcept>
80#include <string>
81#include <vector>
82
83namespace Opm {
84
91template <class TypeTag>
92class FlowProblem : public GetPropType<TypeTag, Properties::BaseProblem>
93 , public FlowGenericProblem<GetPropType<TypeTag, Properties::GridView>,
94 GetPropType<TypeTag, Properties::FluidSystem>>
95{
96protected:
101
110
111 // Grid and world dimension
112 enum { dim = GridView::dimension };
113 enum { dimWorld = GridView::dimensionworld };
114
115 // copy some indices for convenience
116 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
117 enum { numPhases = FluidSystem::numPhases };
118 enum { numComponents = FluidSystem::numComponents };
119
120 enum { enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>() };
121 enum { enableBrine = getPropValue<TypeTag, Properties::EnableBrine>() };
122 enum { enableConvectiveMixing = getPropValue<TypeTag, Properties::EnableConvectiveMixing>() };
123 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
124 enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
125 static constexpr EnergyModules energyModuleType = getPropValue<TypeTag, Properties::EnergyModuleType>();
126 enum { enableFullyImplicitThermal = getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal };
127 enum { enableExperiments = getPropValue<TypeTag, Properties::EnableExperiments>() };
128 enum { enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>() };
129 enum { enableFoam = getPropValue<TypeTag, Properties::EnableFoam>() };
130 enum { enableMICP = Indices::enableMICP };
131 enum { enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>() };
132 enum { enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>() };
133 enum { enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>() };
134 enum { enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>() };
135 enum { enableThermalFluxBoundaries = getPropValue<TypeTag, Properties::EnableThermalFluxBoundaries>() };
136
137 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
138 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
139 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
140
141 // TODO: later, gasCompIdx, oilCompIdx and waterCompIdx should go to the FlowProblemBlackoil in the future
142 // we do not want them in the compositional setting
143 enum { gasCompIdx = FluidSystem::gasCompIdx };
144 enum { oilCompIdx = FluidSystem::oilCompIdx };
145 enum { waterCompIdx = FluidSystem::waterCompIdx };
146
150 using Element = typename GridView::template Codim<0>::Entity;
154 using MaterialLawParams = typename EclMaterialLawManager::MaterialLawParams;
155 using SolidEnergyLawParams = typename EclThermalLawManager::SolidEnergyLawParams;
156 using ThermalConductionLawParams = typename EclThermalLawManager::ThermalConductionLawParams;
163
164 using Toolbox = MathToolbox<Evaluation>;
165 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
166
169 using DirectionalMobilityPtr = Utility::CopyablePtr<DirectionalMobility<TypeTag>>;
170
171public:
177 using BaseType::lame;
180 using BaseType::porosity;
181
185 static void registerParameters()
186 {
187 ParentType::registerParameters();
188
189 registerFlowProblemParameters<Scalar>();
190 }
191
201 static int handlePositionalParameter(std::function<void(const std::string&,
202 const std::string&)> addKey,
203 std::set<std::string>& seenParams,
204 std::string& errorMsg,
205 int,
206 const char** argv,
207 int paramIdx,
208 int)
209 {
210 return detail::eclPositionalParameter(addKey,
211 seenParams,
212 errorMsg,
213 argv,
214 paramIdx);
215 }
216
220 explicit FlowProblem(Simulator& simulator)
221 : ParentType(simulator)
222 , BaseType(simulator.vanguard().eclState(),
223 simulator.vanguard().schedule(),
224 simulator.vanguard().gridView())
225 , transmissibilities_(simulator.vanguard().eclState(),
226 simulator.vanguard().gridView(),
227 simulator.vanguard().cartesianIndexMapper(),
228 simulator.vanguard().grid(),
229 simulator.vanguard().cellCentroids(),
230 (energyModuleType == EnergyModules::FullyImplicitThermal ||
231 energyModuleType == EnergyModules::SequentialImplicitThermal), enableDiffusion,
233 , wellModel_(simulator)
234 , aquiferModel_(simulator)
235 , pffDofData_(simulator.gridView(), this->elementMapper())
236 , tracerModel_(simulator)
237 , temperatureModel_(simulator)
238 {
239 if (! Parameters::Get<Parameters::CheckSatfuncConsistency>()) {
240 // User did not enable the "new" saturation function consistency
241 // check module. Run the original checker instead. This is a
242 // temporary measure.
243 RelpermDiagnostics relpermDiagnostics{};
244 relpermDiagnostics.diagnosis(simulator.vanguard().eclState(),
245 simulator.vanguard().levelCartesianIndexMapper());
246 }
247
248 if (energyModuleType == EnergyModules::SequentialImplicitThermal) {
249 this->enableDriftCompensationTemp_ = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
250 }
251
252 }
253
254 virtual ~FlowProblem() = default;
255
256 void prefetch(const Element& elem) const
257 { this->pffDofData_.prefetch(elem); }
258
270 template <class Restarter>
271 void deserialize(Restarter& res)
272 {
273 // reload the current episode/report step from the deck
274 this->beginEpisode();
275
276 // deserialize the wells
277 wellModel_.deserialize(res);
278
279 // deserialize the aquifer
280 aquiferModel_.deserialize(res);
281 }
282
289 template <class Restarter>
290 void serialize(Restarter& res)
291 {
292 wellModel_.serialize(res);
293
294 aquiferModel_.serialize(res);
295 }
296
297 int episodeIndex() const
298 {
299 return std::max(this->simulator().episodeIndex(), 0);
300 }
301
305 virtual void beginEpisode()
306 {
307 OPM_TIMEBLOCK(beginEpisode);
308 // Proceed to the next report step
309 auto& simulator = this->simulator();
310 int episodeIdx = simulator.episodeIndex();
311 auto& eclState = simulator.vanguard().eclState();
312 const auto& schedule = simulator.vanguard().schedule();
313 const auto& events = schedule[episodeIdx].events();
314
315 if (episodeIdx >= 0 && events.hasEvent(ScheduleEvents::GEO_MODIFIER)) {
316 // bring the contents of the keywords to the current state of the SCHEDULE
317 // section.
318 //
319 // TODO (?): make grid topology changes possible (depending on what exactly
320 // has changed, the grid may need be re-created which has some serious
321 // implications on e.g., the solution of the simulation.)
322 const auto& miniDeck = schedule[episodeIdx].geo_keywords();
323 const auto& cc = simulator.vanguard().grid().comm();
324 eclState.apply_schedule_keywords( miniDeck );
325 eclBroadcast(cc, eclState.getTransMult() );
326
327 // Re-ordering in case of ALUGrid
328 std::function<unsigned int(unsigned int)> equilGridToGrid = [&simulator](unsigned int i) {
329 return simulator.vanguard().gridEquilIdxToGridIdx(i);
330 };
331
332 // re-compute all quantities which may possibly be affected.
333 using TransUpdateQuantities = typename Vanguard::TransmissibilityType::TransUpdateQuantities;
334 transmissibilities_.update(true, TransUpdateQuantities::All, equilGridToGrid);
335 this->referencePorosity_[1] = this->referencePorosity_[0];
337 this->rockFraction_[1] = this->rockFraction_[0];
340 this->model().linearizer().updateDiscretizationParameters();
341 }
342
343 bool tuningEvent = this->beginEpisode_(enableExperiments, this->episodeIndex());
344
345 // set up the wells for the next episode.
346 wellModel_.beginEpisode();
347
348 // set up the aquifers for the next episode.
349 aquiferModel_.beginEpisode();
350
351 // set the size of the initial time step of the episode
352 Scalar dt = limitNextTimeStepSize_(simulator.episodeLength());
353 // negative value of initialTimeStepSize_ indicates no active limit from TSINIT or NEXTSTEP
354 if ( (episodeIdx == 0 || tuningEvent) && this->initialTimeStepSize_ > 0)
355 // allow the size of the initial time step to be set via an external parameter
356 // if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
357 dt = std::min(dt, this->initialTimeStepSize_);
358 simulator.setTimeStepSize(dt);
359 }
360
364 virtual void beginTimeStep()
365 {
366 OPM_TIMEBLOCK(beginTimeStep);
367 const int episodeIdx = this->episodeIndex();
368 const int timeStepSize = this->simulator().timeStepSize();
369
371 episodeIdx,
372 this->simulator().timeStepIndex(),
373 this->simulator().startTime(),
374 this->simulator().time(),
375 timeStepSize,
376 this->simulator().endTime());
377
378 // update maximum water saturation and minimum pressure
379 // used when ROCKCOMP is activated
380 // Do not update max RS first step after a restart
381 this->updateExplicitQuantities_(episodeIdx, timeStepSize, first_step_ && (episodeIdx > 0));
382 first_step_ = false;
383
385 this->model().linearizer().updateBoundaryConditionData();
386 }
387
388 wellModel_.beginTimeStep();
389 aquiferModel_.beginTimeStep();
390 tracerModel_.beginTimeStep();
391 temperatureModel_.beginTimeStep();
392
393 }
394
399 {
400 OPM_TIMEBLOCK(beginIteration);
401 wellModel_.beginIteration();
402 aquiferModel_.beginIteration();
403 }
404
409 {
410 OPM_TIMEBLOCK(endIteration);
411 wellModel_.endIteration();
412 aquiferModel_.endIteration();
413 }
414
418 virtual void endTimeStep()
419 {
420 OPM_TIMEBLOCK(endTimeStep);
421
422#ifndef NDEBUG
423 if constexpr (getPropValue<TypeTag, Properties::EnableDebuggingChecks>()) {
424 // in debug mode, we don't care about performance, so we check
425 // if the model does the right thing (i.e., the mass change
426 // inside the whole reservoir must be equivalent to the fluxes
427 // over the grid's boundaries plus the source rates specified by
428 // the problem).
429 const int rank = this->simulator().gridView().comm().rank();
430 if (rank == 0) {
431 std::cout << "checking conservativeness of solution\n";
432 }
433
434 this->model().checkConservativeness(/*tolerance=*/-1, /*verbose=*/true);
435 if (rank == 0) {
436 std::cout << "solution is sufficiently conservative\n";
437 }
438 }
439#endif // NDEBUG
440
441 auto& simulator = this->simulator();
442 simulator.setTimeStepIndex(simulator.timeStepIndex()+1);
443
444 this->wellModel_.endTimeStep();
445 this->aquiferModel_.endTimeStep();
446 this->tracerModel_.endTimeStep();
447
448 // Compute flux for output
449 this->model().linearizer().updateFlowsInfo();
450
452 OPM_TIMEBLOCK(driftCompansation);
453
454 const auto& residual = this->model().linearizer().residual();
455
456 for (unsigned globalDofIdx = 0; globalDofIdx < residual.size(); globalDofIdx ++) {
457 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(globalDofIdx);
458 this->drift_[sfcdofIdx] = residual[sfcdofIdx] * simulator.timeStepSize();
459
460 if constexpr (getPropValue<TypeTag, Properties::UseVolumetricResidual>()) {
461 this->drift_[sfcdofIdx] *= this->model().dofTotalVolume(sfcdofIdx);
462 }
463 }
464 }
465
466 // Drift compensation needs to be updated before calling the temperature equation
467 if constexpr(energyModuleType == EnergyModules::SequentialImplicitThermal) {
468 this->temperatureModel_.endTimeStep(wellModel_.wellState());
469 }
470 }
471
475 virtual void endEpisode()
476 {
477 const int episodeIdx = this->episodeIndex();
478
479 this->wellModel_.endEpisode();
480 this->aquiferModel_.endEpisode();
481
482 const auto& schedule = this->simulator().vanguard().schedule();
483
484 // End simulation when completed.
485 if (episodeIdx + 1 >= static_cast<int>(schedule.size()) - 1) {
486 this->simulator().setFinished(true);
487 return;
488 }
489
490 // Otherwise, start next episode (report step).
491 this->simulator().startNextEpisode(schedule.stepLength(episodeIdx + 1));
492 }
493
498 virtual void writeOutput(bool verbose)
499 {
500 OPM_TIMEBLOCK(problemWriteOutput);
501
502 if (Parameters::Get<Parameters::EnableWriteAllSolutions>() ||
503 this->episodeWillBeOver())
504 {
505 // Create VTK output as needed.
506 ParentType::writeOutput(verbose);
507 }
508 }
509
513 template <class Context>
514 const DimMatrix& intrinsicPermeability(const Context& context,
515 unsigned spaceIdx,
516 unsigned timeIdx) const
517 {
518 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
519 return transmissibilities_.permeability(globalSpaceIdx);
520 }
521
528 const DimMatrix& intrinsicPermeability(unsigned globalElemIdx) const
529 { return transmissibilities_.permeability(globalElemIdx); }
530
534 template <class Context>
535 Scalar transmissibility(const Context& context,
536 [[maybe_unused]] unsigned fromDofLocalIdx,
537 unsigned toDofLocalIdx) const
538 {
539 assert(fromDofLocalIdx == 0);
540 return pffDofData_.get(context.element(), toDofLocalIdx).transmissibility;
541 }
542
546 Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
547 {
548 return transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
549 }
550
554 template <class Context>
555 Scalar diffusivity(const Context& context,
556 [[maybe_unused]] unsigned fromDofLocalIdx,
557 unsigned toDofLocalIdx) const
558 {
559 assert(fromDofLocalIdx == 0);
560 return *pffDofData_.get(context.element(), toDofLocalIdx).diffusivity;
561 }
562
566 Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
567 return transmissibilities_.diffusivity(globalCellIn, globalCellOut);
568 }
569
573 Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
574 return transmissibilities_.dispersivity(globalCellIn, globalCellOut);
575 }
576
580 Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx,
581 const unsigned boundaryFaceIdx) const
582 {
583 return transmissibilities_.thermalTransmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
584 }
585
586
587
588
592 template <class Context>
593 Scalar transmissibilityBoundary(const Context& elemCtx,
594 unsigned boundaryFaceIdx) const
595 {
596 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
597 return transmissibilities_.transmissibilityBoundary(elemIdx, boundaryFaceIdx);
598 }
599
603 Scalar transmissibilityBoundary(const unsigned globalSpaceIdx,
604 const unsigned boundaryFaceIdx) const
605 {
606 return transmissibilities_.transmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
607 }
608
609
613 Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn,
614 const unsigned globalSpaceIdxOut) const
615 {
616 return transmissibilities_.thermalHalfTrans(globalSpaceIdxIn,globalSpaceIdxOut);
617 }
618
622 template <class Context>
623 Scalar thermalHalfTransmissibilityIn(const Context& context,
624 unsigned faceIdx,
625 unsigned timeIdx) const
626 {
627 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
628 unsigned toDofLocalIdx = face.exteriorIndex();
629 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransIn;
630 }
631
635 template <class Context>
637 unsigned faceIdx,
638 unsigned timeIdx) const
639 {
640 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
641 unsigned toDofLocalIdx = face.exteriorIndex();
642 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransOut;
643 }
644
648 template <class Context>
650 unsigned boundaryFaceIdx) const
651 {
652 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
653 return transmissibilities_.thermalHalfTransBoundary(elemIdx, boundaryFaceIdx);
654 }
655
659 const typename Vanguard::TransmissibilityType& eclTransmissibilities() const
660 { return transmissibilities_; }
661
662
664 { return tracerModel_; }
665
667 { return tracerModel_; }
668
669 TemperatureModel& temperatureModel() // need for restart
670 { return temperatureModel_; }
671
680 template <class Context>
681 Scalar porosity(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
682 {
683 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
684 return this->porosity(globalSpaceIdx, timeIdx);
685 }
686
693 template <class Context>
694 Scalar dofCenterDepth(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
695 {
696 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
697 return this->dofCenterDepth(globalSpaceIdx);
698 }
699
706 Scalar dofCenterDepth(unsigned globalSpaceIdx) const
707 {
708 return this->simulator().vanguard().cellCenterDepth(globalSpaceIdx);
709 }
710
714 template <class Context>
715 Scalar rockCompressibility(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
716 {
717 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
718 return this->rockCompressibility(globalSpaceIdx);
719 }
720
724 template <class Context>
725 Scalar rockBiotComp(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
726 {
727 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
728 return this->rockBiotComp(globalSpaceIdx);
729 }
730
734 template <class Context>
735 Scalar lame(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
736 {
737 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
738 return this->lame(globalSpaceIdx);
739 }
740
744 template <class Context>
745 Scalar biotCoeff(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
746 {
747 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
748 return this->biotCoeff(globalSpaceIdx);
749 }
750
754 template <class Context>
755 Scalar rockReferencePressure(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
756 {
757 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
758 return rockReferencePressure(globalSpaceIdx);
759 }
760
764 Scalar rockReferencePressure(unsigned globalSpaceIdx) const
765 {
766 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
767 if (rock_config.store()) {
768 return asImp_().initialFluidState(globalSpaceIdx).pressure(refPressurePhaseIdx_());
769 }
770 else {
771 if (this->rockParams_.empty())
772 return 1e5;
773
774 unsigned tableIdx = 0;
775 if (!this->rockTableIdx_.empty()) {
776 tableIdx = this->rockTableIdx_[globalSpaceIdx];
777 }
778 return this->rockParams_[tableIdx].referencePressure;
779 }
780 }
781
785 template <class Context>
786 const MaterialLawParams& materialLawParams(const Context& context,
787 unsigned spaceIdx, unsigned timeIdx) const
788 {
789 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
790 return this->materialLawParams(globalSpaceIdx);
791 }
792
793 const MaterialLawParams& materialLawParams(unsigned globalDofIdx) const
794 {
795 return materialLawManager_->materialLawParams(globalDofIdx);
796 }
797
798 const MaterialLawParams& materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
799 {
800 return materialLawManager_->materialLawParams(globalDofIdx, facedir);
801 }
802
806 template <class Context>
808 solidEnergyLawParams(const Context& context,
809 unsigned spaceIdx,
810 unsigned timeIdx) const
811 {
812 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
813 return thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
814 }
815
819 template <class Context>
821 thermalConductionLawParams(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
822 {
823 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
824 return thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
825 }
826
833 std::shared_ptr<const EclMaterialLawManager> materialLawManager() const
834 { return materialLawManager_; }
835
836 template <class FluidState, class ...Args>
838 std::array<Evaluation,numPhases> &mobility,
840 FluidState &fluidState,
841 unsigned globalSpaceIdx) const
842 {
843 using ContainerT = std::array<Evaluation, numPhases>;
844 OPM_TIMEBLOCK_LOCAL(updateRelperms, Subsystem::SatProps);
845 {
846 // calculate relative permeabilities. note that we store the result into the
847 // mobility_ class attribute. the division by the phase viscosity happens later.
848 const auto& materialParams = materialLawParams(globalSpaceIdx);
849 MaterialLaw::template relativePermeabilities<ContainerT, FluidState, Args...>(mobility, materialParams, fluidState);
850 Valgrind::CheckDefined(mobility);
851 }
852 if (materialLawManager_->hasDirectionalRelperms()
853 || materialLawManager_->hasDirectionalImbnum())
854 {
855 using Dir = FaceDir::DirEnum;
856 constexpr int ndim = 3;
857 dirMob = std::make_unique<DirectionalMobility<TypeTag>>();
858 Dir facedirs[ndim] = {Dir::XPlus, Dir::YPlus, Dir::ZPlus};
859 for (int i = 0; i<ndim; i++) {
860 const auto& materialParams = materialLawParams(globalSpaceIdx, facedirs[i]);
861 auto& mob_array = dirMob->getArray(i);
862 MaterialLaw::template relativePermeabilities<ContainerT, FluidState, Args...>(mob_array, materialParams, fluidState);
863 }
864 }
865 }
866
870 std::shared_ptr<EclMaterialLawManager> materialLawManager()
871 { return materialLawManager_; }
872
877 template <class Context>
878 unsigned pvtRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
879 { return pvtRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
880
885 template <class Context>
886 unsigned satnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
887 { return this->satnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
888
893 template <class Context>
894 unsigned miscnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
895 { return this->miscnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
896
901 template <class Context>
902 unsigned plmixnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
903 { return this->plmixnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
904
905 // TODO: polymer related might need to go to the blackoil side
910 template <class Context>
911 Scalar maxPolymerAdsorption(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
912 { return this->maxPolymerAdsorption(context.globalSpaceIndex(spaceIdx, timeIdx)); }
913
917 std::string name() const
918 { return this->simulator().vanguard().caseName(); }
919
923 template <class Context>
924 Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
925 {
926 // use the initial temperature of the DOF if temperature is not a primary
927 // variable
928 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
929 if constexpr (energyModuleType == EnergyModules::SequentialImplicitThermal)
930 return temperatureModel_.temperature(globalDofIdx);
931
932 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
933 }
934
935
936 Scalar temperature(unsigned globalDofIdx, unsigned /*timeIdx*/) const
937 {
938 // use the initial temperature of the DOF if temperature is not a primary
939 // variable
940 if constexpr (energyModuleType == EnergyModules::SequentialImplicitThermal)
941 return temperatureModel_.temperature(globalDofIdx);
942
943 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
944 }
945
947 solidEnergyLawParams(unsigned globalSpaceIdx,
948 unsigned /*timeIdx*/) const
949 {
950 return this->thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
951 }
953 thermalConductionLawParams(unsigned globalSpaceIdx,
954 unsigned /*timeIdx*/)const
955 {
956 return this->thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
957 }
958
968 Scalar maxOilSaturation(unsigned globalDofIdx) const
969 {
970 if (!this->vapparsActive(this->episodeIndex()))
971 return 0.0;
972
973 return this->maxOilSaturation_[globalDofIdx];
974 }
975
985 void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
986 {
987 if (!this->vapparsActive(this->episodeIndex()))
988 return;
989
990 this->maxOilSaturation_[globalDofIdx] = value;
991 }
992
997 {
998 // Calculate all intensive quantities.
999 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/0);
1000
1001 // We also need the intensive quantities for timeIdx == 1
1002 // corresponding to the start of the current timestep, if we
1003 // do not use the storage cache, or if we cannot recycle the
1004 // first iteration storage.
1005 if (!this->model().enableStorageCache() || !this->recycleFirstIterationStorage()) {
1006 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/1);
1007 }
1008
1009 // initialize the wells. Note that this needs to be done after initializing the
1010 // intrinsic permeabilities and the after applying the initial solution because
1011 // the well model uses these...
1012 wellModel_.init();
1013
1014 aquiferModel_.initialSolutionApplied();
1015
1016 const bool invalidateFromHyst = updateHysteresis_();
1017 if (invalidateFromHyst) {
1018 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
1019 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
1020 }
1021 }
1022
1028 template <class Context>
1029 void source(RateVector& rate,
1030 const Context& context,
1031 unsigned spaceIdx,
1032 unsigned timeIdx) const
1033 {
1034 const unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
1035 source(rate, globalDofIdx, timeIdx);
1036 }
1037
1038 void source(RateVector& rate,
1039 unsigned globalDofIdx,
1040 unsigned timeIdx) const
1041 {
1042 OPM_TIMEBLOCK_LOCAL(eclProblemSource, Subsystem::Assembly);
1043 rate = 0.0;
1044
1045 // Add well contribution to source here.
1046 wellModel_.computeTotalRatesForDof(rate, globalDofIdx);
1047
1048 // convert the source term from the total mass rate of the
1049 // cell to the one per unit of volume as used by the model.
1050 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
1051 rate[eqIdx] /= this->model().dofTotalVolume(globalDofIdx);
1052
1053 Valgrind::CheckDefined(rate[eqIdx]);
1054 assert(isfinite(rate[eqIdx]));
1055 }
1056
1057 // Add non-well sources.
1058 addToSourceDense(rate, globalDofIdx, timeIdx);
1059 }
1060
1061 virtual void addToSourceDense(RateVector& rate,
1062 unsigned globalDofIdx,
1063 unsigned timeIdx) const = 0;
1064
1070 const WellModel& wellModel() const
1071 { return wellModel_; }
1072
1074 { return wellModel_; }
1075
1077 { return aquiferModel_; }
1078
1080 { return aquiferModel_; }
1081
1084
1092 {
1093 OPM_TIMEBLOCK(nexTimeStepSize);
1094 // allow external code to do the timestepping
1095 if (this->nextTimeStepSize_ > 0.0)
1096 return this->nextTimeStepSize_;
1097
1098 const auto& simulator = this->simulator();
1099 int episodeIdx = simulator.episodeIndex();
1100
1101 // for the initial episode, we use a fixed time step size
1102 if (episodeIdx < 0)
1103 return this->initialTimeStepSize_;
1104
1105 // ask the newton method for a suggestion. This suggestion will be based on how
1106 // well the previous time step converged. After that, apply the runtime time
1107 // stepping constraints.
1108 const auto& newtonMethod = this->model().newtonMethod();
1109 return limitNextTimeStepSize_(newtonMethod.suggestTimeStepSize(simulator.timeStepSize()));
1110 }
1111
1117 template <class LhsEval>
1118 LhsEval rockCompPoroMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1119 {
1120 OPM_TIMEBLOCK_LOCAL(rockCompPoroMultiplier, Subsystem::PvtProps);
1121 if (this->rockCompPoroMult_.empty() && this->rockCompPoroMultWc_.empty())
1122 return 1.0;
1123
1124 unsigned tableIdx = 0;
1125 if (!this->rockTableIdx_.empty())
1126 tableIdx = this->rockTableIdx_[elementIdx];
1127
1128 const auto& fs = intQuants.fluidState();
1129 LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
1130 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1131 if (!this->minRefPressure_.empty())
1132 // The pore space change is irreversible
1133 effectivePressure =
1134 min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
1135 this->minRefPressure_[elementIdx]);
1136
1137 if (!this->overburdenPressure_.empty())
1138 effectivePressure -= this->overburdenPressure_[elementIdx];
1139
1140 if (rock_config.store()) {
1141 effectivePressure -= asImp_().initialFluidState(elementIdx).pressure(refPressurePhaseIdx_());
1142 }
1143
1144 if (!this->rockCompPoroMult_.empty()) {
1145 return this->rockCompPoroMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1146 }
1147
1148 // water compaction
1149 assert(!this->rockCompPoroMultWc_.empty());
1150 LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1151 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1152
1153 return this->rockCompPoroMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1154 }
1155
1161 template <class LhsEval>
1162 LhsEval rockCompTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1163 {
1164 auto obtain = [](const auto& value)
1165 {
1166 if constexpr (std::is_same_v<LhsEval, Scalar>) {
1167 return getValue(value);
1168 } else {
1169 return value;
1170 }
1171 };
1172 return rockCompTransMultiplier<LhsEval>(intQuants, elementIdx, obtain);
1173 }
1174
1175 template <class LhsEval, class Callback>
1176 LhsEval rockCompTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx, Callback& obtain) const
1177 {
1178 const bool implicit = !this->explicitRockCompaction_;
1179 return implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx, obtain)
1180 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1181 }
1182
1183 template <class LhsEval, class Callback>
1184 LhsEval wellTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx, Callback& obtain) const
1185 {
1186 OPM_TIMEBLOCK_LOCAL(wellTransMultiplier, Subsystem::Wells);
1187
1188 const bool implicit = !this->explicitRockCompaction_;
1189 LhsEval trans_mult = implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx, obtain)
1190 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1191 trans_mult *= this->simulator().problem().template permFactTransMultiplier<LhsEval>(intQuants, elementIdx, obtain);
1192
1193 return trans_mult;
1194 }
1195
1196 std::pair<BCType, RateVector> boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
1197 {
1198 OPM_TIMEBLOCK_LOCAL(boundaryCondition, Subsystem::Assembly);
1200 return { BCType::NONE, RateVector(0.0) };
1201 }
1202 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
1203 const auto& schedule = this->simulator().vanguard().schedule();
1204 if (bcindex_(dir)[globalSpaceIdx] == 0) {
1205 return { BCType::NONE, RateVector(0.0) };
1206 }
1207 if (schedule[this->episodeIndex()].bcprop.size() == 0) {
1208 return { BCType::NONE, RateVector(0.0) };
1209 }
1210 const auto& bc = schedule[this->episodeIndex()].bcprop[bcindex_(dir)[globalSpaceIdx]];
1211 if (bc.bctype!=BCType::RATE) {
1212 return { bc.bctype, RateVector(0.0) };
1213 }
1214
1215 RateVector rate = 0.0;
1216 switch (bc.component) {
1217 case BCComponent::OIL:
1218 rate[FluidSystem::canonicalToActiveCompIdx(oilCompIdx)] = bc.rate;
1219 break;
1220 case BCComponent::GAS:
1221 rate[FluidSystem::canonicalToActiveCompIdx(gasCompIdx)] = bc.rate;
1222 break;
1223 case BCComponent::WATER:
1224 rate[FluidSystem::canonicalToActiveCompIdx(waterCompIdx)] = bc.rate;
1225 break;
1226 case BCComponent::SOLVENT:
1227 this->handleSolventBC(bc, rate);
1228 break;
1229 case BCComponent::POLYMER:
1230 this->handlePolymerBC(bc, rate);
1231 break;
1232 case BCComponent::MICR:
1233 this->handleMicrBC(bc, rate);
1234 break;
1235 case BCComponent::OXYG:
1236 this->handleOxygBC(bc, rate);
1237 break;
1238 case BCComponent::UREA:
1239 this->handleUreaBC(bc, rate);
1240 break;
1241 case BCComponent::NONE:
1242 throw std::logic_error("you need to specify the component when RATE type is set in BC");
1243 break;
1244 }
1245 //TODO add support for enthalpy rate
1246 return {bc.bctype, rate};
1247 }
1248
1249
1250 template<class Serializer>
1251 void serializeOp(Serializer& serializer)
1252 {
1253 serializer(static_cast<BaseType&>(*this));
1254 serializer(drift_);
1255 serializer(wellModel_);
1256 serializer(aquiferModel_);
1257 serializer(tracerModel_);
1258 serializer(*materialLawManager_);
1259 }
1260
1261 const GlobalEqVector& drift() const
1262 {
1263 return drift_;
1264 }
1265
1266private:
1267 Implementation& asImp_()
1268 { return *static_cast<Implementation *>(this); }
1269
1270 const Implementation& asImp_() const
1271 { return *static_cast<const Implementation *>(this); }
1272
1273protected:
1274 template<class UpdateFunc>
1275 void updateProperty_(const std::string& failureMsg,
1276 UpdateFunc func)
1277 {
1278 OPM_TIMEBLOCK(updateProperty);
1279 const auto& model = this->simulator().model();
1280 const auto& primaryVars = model.solution(/*timeIdx*/0);
1281 const auto& vanguard = this->simulator().vanguard();
1282 std::size_t numGridDof = primaryVars.size();
1284#ifdef _OPENMP
1285#pragma omp parallel for
1286#endif
1287 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
1288 const auto& iq = *model.cachedIntensiveQuantities(dofIdx, /*timeIdx=*/ 0);
1289 func(dofIdx, iq);
1290 }
1291 OPM_END_PARALLEL_TRY_CATCH(failureMsg, vanguard.grid().comm());
1292 }
1293
1295 {
1296 OPM_TIMEBLOCK(updateMaxOilSaturation);
1297 int episodeIdx = this->episodeIndex();
1298
1299 // we use VAPPARS
1300 if (this->vapparsActive(episodeIdx)) {
1301 this->updateProperty_("FlowProblem::updateMaxOilSaturation_() failed:",
1302 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1303 {
1304 this->updateMaxOilSaturation_(compressedDofIdx,iq);
1305 });
1306 return true;
1307 }
1308
1309 return false;
1310 }
1311
1312 bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1313 {
1314 OPM_TIMEBLOCK_LOCAL(updateMaxOilSaturation, Subsystem::SatProps);
1315 const auto& fs = iq.fluidState();
1316 const Scalar So = decay<Scalar>(fs.saturation(refPressurePhaseIdx_()));
1317 auto& mos = this->maxOilSaturation_;
1318 if(mos[compressedDofIdx] < So){
1319 mos[compressedDofIdx] = So;
1320 return true;
1321 }else{
1322 return false;
1323 }
1324 }
1325
1327 {
1328 OPM_TIMEBLOCK(updateMaxWaterSaturation);
1329 // water compaction is activated in ROCKCOMP
1330 if (this->maxWaterSaturation_.empty())
1331 return false;
1332
1333 this->maxWaterSaturation_[/*timeIdx=*/1] = this->maxWaterSaturation_[/*timeIdx=*/0];
1334 this->updateProperty_("FlowProblem::updateMaxWaterSaturation_() failed:",
1335 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1336 {
1337 this->updateMaxWaterSaturation_(compressedDofIdx,iq);
1338 });
1339 return true;
1340 }
1341
1342
1343 bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1344 {
1345 OPM_TIMEBLOCK_LOCAL(updateMaxWaterSaturation, Subsystem::SatProps);
1346 const auto& fs = iq.fluidState();
1347 const Scalar Sw = decay<Scalar>(fs.saturation(waterPhaseIdx));
1348 auto& mow = this->maxWaterSaturation_;
1349 if(mow[compressedDofIdx]< Sw){
1350 mow[compressedDofIdx] = Sw;
1351 return true;
1352 }else{
1353 return false;
1354 }
1355 }
1356
1358 {
1359 OPM_TIMEBLOCK(updateMinPressure);
1360 // IRREVERS option is used in ROCKCOMP
1361 if (this->minRefPressure_.empty())
1362 return false;
1363
1364 this->updateProperty_("FlowProblem::updateMinPressure_() failed:",
1365 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1366 {
1367 this->updateMinPressure_(compressedDofIdx,iq);
1368 });
1369 return true;
1370 }
1371
1372 bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities& iq){
1373 OPM_TIMEBLOCK_LOCAL(updateMinPressure, Subsystem::PvtProps);
1374 const auto& fs = iq.fluidState();
1375 const Scalar min_pressure = getValue(fs.pressure(refPressurePhaseIdx_()));
1376 auto& min_pressures = this->minRefPressure_;
1377 if(min_pressures[compressedDofIdx]> min_pressure){
1378 min_pressures[compressedDofIdx] = min_pressure;
1379 return true;
1380 }else{
1381 return false;
1382 }
1383 }
1384
1385 // \brief Function to assign field properties of type double, on the leaf grid view.
1386 //
1387 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1388 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1389 std::function<std::vector<double>(const FieldPropsManager&, const std::string&)>
1391 {
1392 const auto& lookup = this->lookUpData_;
1393 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString)
1394 {
1395 return lookup.assignFieldPropsDoubleOnLeaf(fieldPropManager, propString);
1396 };
1397 }
1398
1399 // \brief Function to assign field properties of type int, unsigned int, ..., on the leaf grid view.
1400 //
1401 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1402 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1403 template<typename IntType>
1404 std::function<std::vector<IntType>(const FieldPropsManager&, const std::string&, bool)>
1406 {
1407 const auto& lookup = this->lookUpData_;
1408 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString, bool needsTranslation)
1409 {
1410 return lookup.template assignFieldPropsIntOnLeaf<IntType>(fieldPropManager, propString, needsTranslation);
1411 };
1412 }
1413
1415 {
1416 OPM_TIMEBLOCK(readMaterialParameters);
1417 const auto& simulator = this->simulator();
1418 const auto& vanguard = simulator.vanguard();
1419 const auto& eclState = vanguard.eclState();
1420
1421 // the PVT and saturation region numbers
1423 this->updatePvtnum_();
1424 this->updateSatnum_();
1425
1426 // the MISC region numbers (solvent model)
1427 this->updateMiscnum_();
1428 // the PLMIX region numbers (polymer model)
1429 this->updatePlmixnum_();
1430
1431 OPM_END_PARALLEL_TRY_CATCH("Invalid region numbers: ", vanguard.gridView().comm());
1433 // porosity
1435 this->referencePorosity_[1] = this->referencePorosity_[0];
1437
1439 // rock fraction
1441 this->rockFraction_[1] = this->rockFraction_[0];
1443
1445 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1446 materialLawManager_ = std::make_shared<EclMaterialLawManager>();
1447 materialLawManager_->initFromState(eclState);
1448 materialLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1449 this-> template fieldPropIntTypeOnLeafAssigner_<int>(),
1452 }
1453
1455 {
1456 if constexpr (energyModuleType == EnergyModules::FullyImplicitThermal ||
1457 energyModuleType == EnergyModules::SequentialImplicitThermal )
1458 {
1459 const auto& simulator = this->simulator();
1460 const auto& vanguard = simulator.vanguard();
1461 const auto& eclState = vanguard.eclState();
1462
1463 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1464 thermalLawManager_ = std::make_shared<EclThermalLawManager>();
1465 thermalLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1467 this-> template fieldPropIntTypeOnLeafAssigner_<unsigned int>());
1468 }
1469 }
1470
1472 {
1473 const auto& simulator = this->simulator();
1474 const auto& vanguard = simulator.vanguard();
1475 const auto& eclState = vanguard.eclState();
1476
1477 std::size_t numDof = this->model().numGridDof();
1478
1479 this->referencePorosity_[/*timeIdx=*/0].resize(numDof);
1480
1481 const auto& fp = eclState.fieldProps();
1482 const std::vector<double> porvData = this -> fieldPropDoubleOnLeafAssigner_()(fp, "PORV");
1483 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1484 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(dofIdx);
1485 Scalar poreVolume = porvData[dofIdx];
1486
1487 // we define the porosity as the accumulated pore volume divided by the
1488 // geometric volume of the element. Note that -- in pathetic cases -- it can
1489 // be larger than 1.0!
1490 Scalar dofVolume = simulator.model().dofTotalVolume(sfcdofIdx);
1491 assert(dofVolume > 0.0);
1492 this->referencePorosity_[/*timeIdx=*/0][sfcdofIdx] = poreVolume/dofVolume;
1493 }
1494 }
1495
1497
1498 const bool solveEnergyEquation = (energyModuleType == EnergyModules::FullyImplicitThermal ||
1499 energyModuleType == EnergyModules::SequentialImplicitThermal);
1500 if (!solveEnergyEquation)
1501 return;
1502
1503 const auto& simulator = this->simulator();
1504 const auto& vanguard = simulator.vanguard();
1505 const auto& eclState = vanguard.eclState();
1506
1507 std::size_t numDof = this->model().numGridDof();
1508 this->rockFraction_[/*timeIdx=*/0].resize(numDof);
1509 // For the energy equation, we need the volume of the rock.
1510 // The volume of the rock is computed by rockFraction * geometric volume of the element.
1511 // The reference porosity is defined as porosity * ntg * pore-volume-multiplier.
1512 // A common practice in reservoir simulation is to use large pore-volume-multipliers in boundary cells
1513 // to model boundary conditions other than no-flow. This may result in reference porosities that are larger than 1.
1514 // A simple (1-reference porosity) * geometric volume of the element may give unphysical results.
1515 // We therefore instead consider the pore-volume-multiplier as a volume multiplier. The rock fraction is thus given by
1516 // (1 - porosity * ntg) * pore-volume-multiplier = (1 - porosity * ntg) * reference porosity / (porosity * ntg)
1517 const auto& fp = eclState.fieldProps();
1518 const std::vector<double> poroData = this->fieldPropDoubleOnLeafAssigner_()(fp, "PORO");
1519 const std::vector<double> ntgData = this->fieldPropDoubleOnLeafAssigner_()(fp, "NTG");
1520
1521 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1522 const auto ntg = ntgData[dofIdx];
1523 const auto poro_eff = ntg * poroData[dofIdx];
1524 const int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(dofIdx);
1525 const auto rock_fraction = (1 - poro_eff) * this->referencePorosity_[/*timeIdx=*/0][sfcdofIdx] / poro_eff;
1526 this->rockFraction_[/*timeIdx=*/0][sfcdofIdx] = rock_fraction;
1527 }
1528 }
1529
1531 {
1532 // TODO: whether we should move this to FlowProblemBlackoil
1533 const auto& simulator = this->simulator();
1534 const auto& vanguard = simulator.vanguard();
1535 const auto& eclState = vanguard.eclState();
1536
1537 if (eclState.getInitConfig().hasEquil())
1539 else
1541
1542 //initialize min/max values
1543 std::size_t numElems = this->model().numGridDof();
1544 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1545 const auto& fs = asImp_().initialFluidStates()[elemIdx];
1546 if (!this->maxWaterSaturation_.empty() && waterPhaseIdx > -1)
1547 this->maxWaterSaturation_[elemIdx] = std::max(this->maxWaterSaturation_[elemIdx], fs.saturation(waterPhaseIdx));
1548 if (!this->maxOilSaturation_.empty() && oilPhaseIdx > -1)
1549 this->maxOilSaturation_[elemIdx] = std::max(this->maxOilSaturation_[elemIdx], fs.saturation(oilPhaseIdx));
1550 if (!this->minRefPressure_.empty() && refPressurePhaseIdx_() > -1)
1551 this->minRefPressure_[elemIdx] = std::min(this->minRefPressure_[elemIdx], fs.pressure(refPressurePhaseIdx_()));
1552 }
1553 }
1554
1555 virtual void readEquilInitialCondition_() = 0;
1557
1558 // update the hysteresis parameters of the material laws for the whole grid
1560 {
1561 if (!materialLawManager_->enableHysteresis())
1562 return false;
1563
1564 // we need to update the hysteresis data for _all_ elements (i.e., not just the
1565 // interior ones) to avoid desynchronization of the processes in the parallel case!
1566 this->updateProperty_("FlowProblem::updateHysteresis_() failed:",
1567 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1568 {
1569 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1570 });
1571 return true;
1572 }
1573
1574
1575 bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1576 {
1577 OPM_TIMEBLOCK_LOCAL(updateHysteresis_, Subsystem::SatProps);
1578 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1579 //TODO change materials to give a bool
1580 return true;
1581 }
1582
1583 Scalar getRockCompTransMultVal(std::size_t dofIdx) const
1584 {
1585 if (this->rockCompTransMultVal_.empty())
1586 return 1.0;
1587
1588 return this->rockCompTransMultVal_[dofIdx];
1589 }
1590
1591protected:
1593 {
1594 ConditionalStorage<enableFullyImplicitThermal, Scalar> thermalHalfTransIn;
1595 ConditionalStorage<enableFullyImplicitThermal, Scalar> thermalHalfTransOut;
1596 ConditionalStorage<enableDiffusion, Scalar> diffusivity;
1597 ConditionalStorage<enableDispersion, Scalar> dispersivity;
1599 };
1600
1601 // update the prefetch friendly data object
1603 {
1604 const auto& distFn =
1605 [this](PffDofData_& dofData,
1606 const Stencil& stencil,
1607 unsigned localDofIdx)
1608 -> void
1609 {
1610 const auto& elementMapper = this->model().elementMapper();
1611
1612 unsigned globalElemIdx = elementMapper.index(stencil.entity(localDofIdx));
1613 if (localDofIdx != 0) {
1614 unsigned globalCenterElemIdx = elementMapper.index(stencil.entity(/*dofIdx=*/0));
1615 dofData.transmissibility = transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
1616
1617 if constexpr (enableFullyImplicitThermal) {
1618 *dofData.thermalHalfTransIn = transmissibilities_.thermalHalfTrans(globalCenterElemIdx, globalElemIdx);
1619 *dofData.thermalHalfTransOut = transmissibilities_.thermalHalfTrans(globalElemIdx, globalCenterElemIdx);
1620 }
1621 if constexpr (enableDiffusion)
1622 *dofData.diffusivity = transmissibilities_.diffusivity(globalCenterElemIdx, globalElemIdx);
1623 if (enableDispersion)
1624 dofData.dispersivity = transmissibilities_.dispersivity(globalCenterElemIdx, globalElemIdx);
1625 }
1626 };
1627
1628 pffDofData_.update(distFn);
1629 }
1630
1631 virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart) = 0;
1632
1634 {
1635 const auto& simulator = this->simulator();
1636 const auto& vanguard = simulator.vanguard();
1637 const auto& bcconfig = vanguard.eclState().getSimulationConfig().bcconfig();
1638 if (bcconfig.size() > 0) {
1640
1641 std::size_t numCartDof = vanguard.cartesianSize();
1642 unsigned numElems = vanguard.gridView().size(/*codim=*/0);
1643 std::vector<int> cartesianToCompressedElemIdx(numCartDof, -1);
1644
1645 for (unsigned elemIdx = 0; elemIdx < numElems; ++elemIdx)
1646 cartesianToCompressedElemIdx[vanguard.cartesianIndex(elemIdx)] = elemIdx;
1647
1648 bcindex_.resize(numElems, 0);
1649 auto loopAndApply = [&cartesianToCompressedElemIdx,
1650 &vanguard](const auto& bcface,
1651 auto apply)
1652 {
1653 for (int i = bcface.i1; i <= bcface.i2; ++i) {
1654 for (int j = bcface.j1; j <= bcface.j2; ++j) {
1655 for (int k = bcface.k1; k <= bcface.k2; ++k) {
1656 std::array<int, 3> tmp = {i,j,k};
1657 auto elemIdx = cartesianToCompressedElemIdx[vanguard.cartesianIndex(tmp)];
1658 if (elemIdx >= 0)
1659 apply(elemIdx);
1660 }
1661 }
1662 }
1663 };
1664 for (const auto& bcface : bcconfig) {
1665 std::vector<int>& data = bcindex_(bcface.dir);
1666 const int index = bcface.index;
1667 loopAndApply(bcface,
1668 [&data,index](int elemIdx)
1669 { data[elemIdx] = index; });
1670 }
1671 }
1672 }
1673
1674 // this method applies the runtime constraints specified via the deck and/or command
1675 // line parameters for the size of the next time step.
1677 {
1678 if constexpr (enableExperiments) {
1679 const auto& simulator = this->simulator();
1680 const auto& schedule = simulator.vanguard().schedule();
1681 int episodeIdx = simulator.episodeIndex();
1682
1683 // first thing in the morning, limit the time step size to the maximum size
1684 Scalar maxTimeStepSize = Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60;
1685 int reportStepIdx = std::max(episodeIdx, 0);
1686 if (this->enableTuning_) {
1687 const auto& tuning = schedule[reportStepIdx].tuning();
1688 maxTimeStepSize = tuning.TSMAXZ;
1689 }
1690
1691 dtNext = std::min(dtNext, maxTimeStepSize);
1692
1693 Scalar remainingEpisodeTime =
1694 simulator.episodeStartTime() + simulator.episodeLength()
1695 - (simulator.startTime() + simulator.time());
1696 assert(remainingEpisodeTime >= 0.0);
1697
1698 // if we would have a small amount of time left over in the current episode, make
1699 // two equal time steps instead of a big and a small one
1700 if (remainingEpisodeTime/2.0 < dtNext && dtNext < remainingEpisodeTime*(1.0 - 1e-5))
1701 // note: limiting to the maximum time step size here is probably not strictly
1702 // necessary, but it should not hurt and is more fool-proof
1703 dtNext = std::min(maxTimeStepSize, remainingEpisodeTime/2.0);
1704
1705 if (simulator.episodeStarts()) {
1706 // if a well event occurred, respect the limit for the maximum time step after
1707 // that, too
1708 const auto& events = simulator.vanguard().schedule()[reportStepIdx].events();
1709 bool wellEventOccured =
1710 events.hasEvent(ScheduleEvents::NEW_WELL)
1711 || events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE)
1712 || events.hasEvent(ScheduleEvents::INJECTION_UPDATE)
1713 || events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
1714 if (episodeIdx >= 0 && wellEventOccured && this->maxTimeStepAfterWellEvent_ > 0)
1715 dtNext = std::min(dtNext, this->maxTimeStepAfterWellEvent_);
1716 }
1717 }
1718
1719 return dtNext;
1720 }
1721
1723 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1724 return oilPhaseIdx;
1725 }
1726 else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
1727 return gasPhaseIdx;
1728 }
1729 else {
1730 return waterPhaseIdx;
1731 }
1732 }
1733
1735 {
1736 const auto& model = this->simulator().model();
1737 std::size_t numGridDof = this->model().numGridDof();
1738 this->rockCompTransMultVal_.resize(numGridDof, 1.0);
1739 for (std::size_t elementIdx = 0; elementIdx < numGridDof; ++elementIdx) {
1740 const auto& iq = *model.cachedIntensiveQuantities(elementIdx, /*timeIdx=*/ 0);
1741 Scalar trans_mult = computeRockCompTransMultiplier_<Scalar>(iq, elementIdx);
1742 this->rockCompTransMultVal_[elementIdx] = trans_mult;
1743 }
1744 }
1745
1751 template <class LhsEval>
1752 LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1753 {
1754 auto obtain = [](const auto& value)
1755 {
1756 if constexpr (std::is_same_v<LhsEval, Scalar>) {
1757 return getValue(value);
1758 } else {
1759 return value;
1760 }
1761 };
1762
1763 return computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx, obtain);
1764 }
1765
1766 template <class LhsEval, class Callback>
1767 LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities& intQuants, unsigned elementIdx, Callback& obtain) const
1768 {
1769 OPM_TIMEBLOCK_LOCAL(computeRockCompTransMultiplier, Subsystem::PvtProps);
1770 if (this->rockCompTransMult_.empty() && this->rockCompTransMultWc_.empty())
1771 return 1.0;
1772
1773 unsigned tableIdx = 0;
1774 if (!this->rockTableIdx_.empty())
1775 tableIdx = this->rockTableIdx_[elementIdx];
1776
1777 const auto& fs = intQuants.fluidState();
1778 LhsEval effectivePressure = obtain(fs.pressure(refPressurePhaseIdx_()));
1779 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1780 if (!this->minRefPressure_.empty())
1781 // The pore space change is irreversible
1782 effectivePressure =
1783 min(obtain(fs.pressure(refPressurePhaseIdx_())),
1784 this->minRefPressure_[elementIdx]);
1785
1786 if (!this->overburdenPressure_.empty())
1787 effectivePressure -= this->overburdenPressure_[elementIdx];
1788
1789 if (rock_config.store()) {
1790 effectivePressure -= asImp_().initialFluidState(elementIdx).pressure(refPressurePhaseIdx_());
1791 }
1792
1793 if (!this->rockCompTransMult_.empty())
1794 return this->rockCompTransMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1795
1796 // water compaction
1797 assert(!this->rockCompTransMultWc_.empty());
1798 LhsEval SwMax = max(obtain(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1799 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1800
1801 return this->rockCompTransMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1802 }
1803
1804 typename Vanguard::TransmissibilityType transmissibilities_;
1805
1806 std::shared_ptr<EclMaterialLawManager> materialLawManager_;
1807 std::shared_ptr<EclThermalLawManager> thermalLawManager_;
1808
1810
1813
1817
1818 template<class T>
1819 struct BCData
1820 {
1821 std::array<std::vector<T>,6> data;
1822
1823 void resize(std::size_t size, T defVal)
1824 {
1825 for (auto& d : data)
1826 d.resize(size, defVal);
1827 }
1828
1829 const std::vector<T>& operator()(FaceDir::DirEnum dir) const
1830 {
1831 if (dir == FaceDir::DirEnum::Unknown)
1832 throw std::runtime_error("Tried to access BC data for the 'Unknown' direction");
1833 int idx = 0;
1834 int div = static_cast<int>(dir);
1835 while ((div /= 2) >= 1)
1836 ++idx;
1837 assert(idx >= 0 && idx <= 5);
1838 return data[idx];
1839 }
1840
1841 std::vector<T>& operator()(FaceDir::DirEnum dir)
1842 {
1843 return const_cast<std::vector<T>&>(std::as_const(*this)(dir));
1844 }
1845 };
1846
1847 virtual void handleSolventBC(const BCProp::BCFace&, RateVector&) const = 0;
1848
1849 virtual void handlePolymerBC(const BCProp::BCFace&, RateVector&) const = 0;
1850
1851 virtual void handleMicrBC(const BCProp::BCFace&, RateVector&) const = 0;
1852
1853 virtual void handleOxygBC(const BCProp::BCFace&, RateVector&) const = 0;
1854
1855 virtual void handleUreaBC(const BCProp::BCFace&, RateVector&) const = 0;
1856
1859 bool first_step_ = true;
1860
1863 virtual bool episodeWillBeOver() const
1864 {
1865 return this->simulator().episodeWillBeOver();
1866 }
1867};
1868
1869} // namespace Opm
1870
1871#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:797
unsigned pvtRegionIndex(unsigned elemIdx) const
Returns the index the relevant PVT region given a cell index.
Definition: FlowGenericProblem_impl.hpp:756
std::array< std::vector< Scalar >, 2 > rockFraction_
Definition: FlowGenericProblem.hpp:338
std::function< unsigned(unsigned)> lookupIdxOnLevelZeroAssigner_()
Definition: FlowGenericProblem.hpp:384
std::vector< TabulatedTwoDFunction > rockCompPoroMultWc_
Definition: FlowGenericProblem.hpp:347
Scalar porosity(unsigned globalSpaceIdx, unsigned timeIdx) const
Direct indexed access to the porosity.
Definition: FlowGenericProblem_impl.hpp:338
Scalar rockCompressibility(unsigned globalSpaceIdx) const
Definition: FlowGenericProblem_impl.hpp:323
unsigned miscnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant MISC region given a cell index.
Definition: FlowGenericProblem_impl.hpp:776
unsigned satnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant saturation function region given a cell index.
Definition: FlowGenericProblem_impl.hpp:766
Scalar lame(unsigned elementIdx) const
Direct access to Lame's second parameter in an element.
Definition: FlowGenericProblem_impl.hpp:357
void beginTimeStep_(bool enableExperiments, int episodeIdx, int timeStepIndex, Scalar startTime, Scalar time, Scalar timeStepSize, Scalar endTime)
Definition: FlowGenericProblem_impl.hpp:508
std::vector< TabulatedTwoDFunction > rockCompTransMultWc_
Definition: FlowGenericProblem.hpp:348
unsigned plmixnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant PLMIXNUM (for polymer module) region given a cell index.
Definition: FlowGenericProblem_impl.hpp:786
std::array< std::vector< Scalar >, 2 > referencePorosity_
Definition: FlowGenericProblem.hpp:337
Scalar rockBiotComp(unsigned elementIdx) const
Returns the rock compressibility of an element due to poroelasticity.
Definition: FlowGenericProblem_impl.hpp:346
bool beginEpisode_(bool enableExperiments, int episodeIdx)
Definition: FlowGenericProblem_impl.hpp:471
bool shouldWriteOutput() const
Always returns true. The ecl output writer takes care of the rest.
Definition: FlowGenericProblem.hpp:294
static std::string helpPreamble(int, const char **argv)
Definition: FlowGenericProblem_impl.hpp:115
bool shouldWriteRestartFile() const
Returns true if an eWoms restart file should be written to disk.
Definition: FlowGenericProblem.hpp:303
Scalar biotCoeff(unsigned elementIdx) const
Direct access to Biot coefficient in an element.
Definition: FlowGenericProblem_impl.hpp:388
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowProblem.hpp:95
virtual bool episodeWillBeOver() const
Definition: FlowProblem.hpp:1863
const WellModel & wellModel() const
Returns a reference to the ECL well manager used by the problem.
Definition: FlowProblem.hpp:1070
GetPropType< TypeTag, Properties::Evaluation > Evaluation
Definition: FlowProblem.hpp:159
Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
Direct access to the transmissibility between two elements.
Definition: FlowProblem.hpp:546
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:201
bool nonTrivialBoundaryConditions() const
Definition: FlowProblem.hpp:1082
const GlobalEqVector & drift() const
Definition: FlowProblem.hpp:1261
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:498
std::function< std::vector< IntType >(const FieldPropsManager &, const std::string &, bool)> fieldPropIntTypeOnLeafAssigner_()
Definition: FlowProblem.hpp:1405
const DimMatrix & intrinsicPermeability(unsigned globalElemIdx) const
This method returns the intrinsic permeability tensor given a global element index.
Definition: FlowProblem.hpp:528
@ gasPhaseIdx
Definition: FlowProblem.hpp:137
@ enableFullyImplicitThermal
Definition: FlowProblem.hpp:126
@ enableConvectiveMixing
Definition: FlowProblem.hpp:122
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:878
LhsEval wellTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Definition: FlowProblem.hpp:1184
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:681
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition: FlowProblem.hpp:398
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:108
GetPropType< TypeTag, Properties::DofMapper > DofMapper
Definition: FlowProblem.hpp:158
@ enableBrine
Definition: FlowProblem.hpp:121
Scalar rockBiotComp(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:725
@ enableSolvent
Definition: FlowProblem.hpp:134
@ dim
Definition: FlowProblem.hpp:112
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: FlowProblem.hpp:102
typename EclThermalLawManager::SolidEnergyLawParams SolidEnergyLawParams
Definition: FlowProblem.hpp:155
bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1575
GetPropType< TypeTag, Properties::BaseProblem > ParentType
Definition: FlowProblem.hpp:99
GetPropType< TypeTag, Properties::EqVector > EqVector
Definition: FlowProblem.hpp:107
unsigned satnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:886
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:636
@ numEq
Definition: FlowProblem.hpp:116
bool first_step_
Definition: FlowProblem.hpp:1859
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:151
const ThermalConductionLawParams & thermalConductionLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:953
AquiferModel aquiferModel_
Definition: FlowProblem.hpp:1812
GlobalEqVector drift_
Definition: FlowProblem.hpp:1809
bool updateMinPressure_()
Definition: FlowProblem.hpp:1357
std::function< std::vector< double >(const FieldPropsManager &, const std::string &)> fieldPropDoubleOnLeafAssigner_()
Definition: FlowProblem.hpp:1390
Scalar transmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:593
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:148
void updateReferencePorosity_()
Definition: FlowProblem.hpp:1471
Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn, const unsigned globalSpaceIdxOut) const
Definition: FlowProblem.hpp:613
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Definition: FlowProblem.hpp:1767
BCData< int > bcindex_
Definition: FlowProblem.hpp:1857
GetPropType< TypeTag, Properties::TracerModel > TracerModel
Definition: FlowProblem.hpp:168
@ enableExtbo
Definition: FlowProblem.hpp:128
Scalar rockCompressibility(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:715
bool updateMaxWaterSaturation_()
Definition: FlowProblem.hpp:1326
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: FlowProblem.hpp:165
Scalar maxOilSaturation(unsigned globalDofIdx) const
Returns an element's historic maximum oil phase saturation that was observed during the simulation.
Definition: FlowProblem.hpp:968
bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1372
std::string name() const
The problem name.
Definition: FlowProblem.hpp:917
int episodeIndex() const
Definition: FlowProblem.hpp:297
LhsEval rockCompTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1162
void endIteration()
Called by the simulator after each Newton-Raphson iteration.
Definition: FlowProblem.hpp:408
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:109
FlowProblem(Simulator &simulator)
Definition: FlowProblem.hpp:220
@ enableExperiments
Definition: FlowProblem.hpp:127
@ enableFoam
Definition: FlowProblem.hpp:129
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:133
virtual void handleSolventBC(const BCProp::BCFace &, RateVector &) const =0
@ numPhases
Definition: FlowProblem.hpp:117
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: FlowProblem.hpp:106
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: FlowProblem.hpp:149
const Vanguard::TransmissibilityType & eclTransmissibilities() const
Return a reference to the object that handles the "raw" transmissibilities.
Definition: FlowProblem.hpp:659
void source(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:1038
virtual void readEquilInitialCondition_()=0
Scalar nextTimeStepSize() const
Propose the size of the next time step to the simulator.
Definition: FlowProblem.hpp:1091
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1752
std::pair< BCType, RateVector > boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
Definition: FlowProblem.hpp:1196
@ enableDispersion
Definition: FlowProblem.hpp:124
typename EclMaterialLawManager::MaterialLawParams MaterialLawParams
Definition: FlowProblem.hpp:154
virtual ~FlowProblem()=default
PffGridVector< GridView, Stencil, PffDofData_, DofMapper > pffDofData_
Definition: FlowProblem.hpp:1814
@ enableDiffusion
Definition: FlowProblem.hpp:123
const ThermalConductionLawParams & thermalConductionLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:821
Scalar transmissibility(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:535
TracerModel tracerModel_
Definition: FlowProblem.hpp:1815
virtual void handlePolymerBC(const BCProp::BCFace &, RateVector &) const =0
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:135
const TracerModel & tracerModel() const
Definition: FlowProblem.hpp:663
WellModel wellModel_
Definition: FlowProblem.hpp:1811
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:305
const SolidEnergyLawParams & solidEnergyLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:947
Scalar getRockCompTransMultVal(std::size_t dofIdx) const
Definition: FlowProblem.hpp:1583
virtual void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:364
Scalar dofCenterDepth(unsigned globalSpaceIdx) const
Direct indexed acces to the depth of an degree of freedom [m].
Definition: FlowProblem.hpp:706
typename GetProp< TypeTag, Properties::MaterialLaw >::EclMaterialLawManager EclMaterialLawManager
Definition: FlowProblem.hpp:152
const MaterialLawParams & materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
Definition: FlowProblem.hpp:798
std::shared_ptr< const EclMaterialLawManager > materialLawManager() const
Returns the ECL material law manager.
Definition: FlowProblem.hpp:833
unsigned plmixnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:902
const MaterialLawParams & materialLawParams(unsigned globalDofIdx) const
Definition: FlowProblem.hpp:793
Scalar temperature(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:924
Scalar thermalHalfTransmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:649
unsigned miscnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:894
void updateRockCompTransMultVal_()
Definition: FlowProblem.hpp:1734
std::shared_ptr< EclThermalLawManager > thermalLawManager_
Definition: FlowProblem.hpp:1807
Scalar limitNextTimeStepSize_(Scalar dtNext) const
Definition: FlowProblem.hpp:1676
GetPropType< TypeTag, Properties::Stencil > Stencil
Definition: FlowProblem.hpp:104
typename GetProp< TypeTag, Properties::SolidEnergyLaw >::EclThermalLawManager EclThermalLawManager
Definition: FlowProblem.hpp:153
virtual void readExplicitInitialCondition_()=0
virtual void handleMicrBC(const BCProp::BCFace &, RateVector &) const =0
GetPropType< TypeTag, Properties::WellModel > WellModel
Definition: FlowProblem.hpp:161
const DimMatrix & intrinsicPermeability(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:514
bool updateHysteresis_()
Definition: FlowProblem.hpp:1559
void readThermalParameters_()
Definition: FlowProblem.hpp:1454
@ enableBioeffects
Definition: FlowProblem.hpp:120
void serializeOp(Serializer &serializer)
Definition: FlowProblem.hpp:1251
Scalar maxPolymerAdsorption(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the max polymer adsorption value.
Definition: FlowProblem.hpp:911
@ gasCompIdx
Definition: FlowProblem.hpp:143
Scalar thermalHalfTransmissibilityIn(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:623
void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
Sets an element's maximum oil phase saturation observed during the simulation.
Definition: FlowProblem.hpp:985
@ enablePolymer
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:808
typename EclThermalLawManager::ThermalConductionLawParams ThermalConductionLawParams
Definition: FlowProblem.hpp:156
Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:573
AquiferModel & mutableAquiferModel()
Definition: FlowProblem.hpp:1079
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:160
Scalar temperature(unsigned globalDofIdx, unsigned) const
Definition: FlowProblem.hpp:936
Scalar lame(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:735
GetPropType< TypeTag, Properties::AquiferModel > AquiferModel
Definition: FlowProblem.hpp:162
@ oilPhaseIdx
Definition: FlowProblem.hpp:138
TemperatureModel temperatureModel_
Definition: FlowProblem.hpp:1816
std::shared_ptr< EclMaterialLawManager > materialLawManager_
Definition: FlowProblem.hpp:1806
WellModel & wellModel()
Definition: FlowProblem.hpp:1073
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowProblem.hpp:103
GetPropType< TypeTag, Properties::TemperatureModel > TemperatureModel
Definition: FlowProblem.hpp:167
void updateProperty_(const std::string &failureMsg, UpdateFunc func)
Definition: FlowProblem.hpp:1275
LhsEval rockCompTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Definition: FlowProblem.hpp:1176
bool updateMaxOilSaturation_()
Definition: FlowProblem.hpp:1294
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblem.hpp:185
void updatePffDofData_()
Definition: FlowProblem.hpp:1602
virtual void endEpisode()
Called by the simulator after the end of an episode.
Definition: FlowProblem.hpp:475
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:147
bool nonTrivialBoundaryConditions_
Definition: FlowProblem.hpp:1858
@ numComponents
Definition: FlowProblem.hpp:118
@ enablePolymerMolarWeight
Definition: FlowProblem.hpp:132
GetPropType< TypeTag, Properties::Problem > Implementation
Definition: FlowProblem.hpp:100
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1633
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:837
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1804
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:1029
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:418
TemperatureModel & temperatureModel()
Definition: FlowProblem.hpp:669
Utility::CopyablePtr< DirectionalMobility< TypeTag > > DirectionalMobilityPtr
Definition: FlowProblem.hpp:169
virtual void readInitialCondition_()
Definition: FlowProblem.hpp:1530
virtual void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblem.hpp:996
const MaterialLawParams & materialLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:786
Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:580
@ enableMICP
Definition: FlowProblem.hpp:130
static constexpr EnergyModules energyModuleType
Definition: FlowProblem.hpp:125
Scalar diffusivity(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:555
Scalar rockReferencePressure(unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:764
@ waterPhaseIdx
Definition: FlowProblem.hpp:139
void deserialize(Restarter &res)
This method restores the complete state of the problem and its sub-objects from disk.
Definition: FlowProblem.hpp:271
Scalar biotCoeff(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:745
std::shared_ptr< EclMaterialLawManager > materialLawManager()
Definition: FlowProblem.hpp:870
@ dimWorld
Definition: FlowProblem.hpp:113
bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1343
typename GridView::template Codim< 0 >::Entity Element
Definition: FlowProblem.hpp:150
Scalar dofCenterDepth(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the depth of an degree of freedom [m].
Definition: FlowProblem.hpp:694
@ oilCompIdx
Definition: FlowProblem.hpp:144
bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1312
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: FlowProblem.hpp:105
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: FlowProblem.hpp:157
const AquiferModel & aquiferModel() const
Definition: FlowProblem.hpp:1076
@ waterCompIdx
Definition: FlowProblem.hpp:145
int refPressurePhaseIdx_() const
Definition: FlowProblem.hpp:1722
LhsEval rockCompPoroMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the porosity multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1118
Scalar transmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:603
TracerModel & tracerModel()
Definition: FlowProblem.hpp:666
void readMaterialParameters_()
Definition: FlowProblem.hpp:1414
void serialize(Restarter &res)
This method writes the complete state of the problem and its subobjects to disk.
Definition: FlowProblem.hpp:290
MathToolbox< Evaluation > Toolbox
Definition: FlowProblem.hpp:164
Scalar rockReferencePressure(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:755
Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:566
void updateRockFraction_()
Definition: FlowProblem.hpp:1496
void prefetch(const Element &elem) const
Definition: FlowProblem.hpp:256
A random-access container which stores data attached to a grid's degrees of freedom in a prefetch fri...
Definition: pffgridvector.hh: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: blackoilbioeffectsmodules.hh:45
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:1820
const std::vector< T > & operator()(FaceDir::DirEnum dir) const
Definition: FlowProblem.hpp:1829
void resize(std::size_t size, T defVal)
Definition: FlowProblem.hpp:1823
std::vector< T > & operator()(FaceDir::DirEnum dir)
Definition: FlowProblem.hpp:1841
std::array< std::vector< T >, 6 > data
Definition: FlowProblem.hpp:1821
Definition: FlowProblem.hpp:1593
ConditionalStorage< enableFullyImplicitThermal, Scalar > thermalHalfTransOut
Definition: FlowProblem.hpp:1595
ConditionalStorage< enableFullyImplicitThermal, Scalar > thermalHalfTransIn
Definition: FlowProblem.hpp:1594
ConditionalStorage< enableDiffusion, Scalar > diffusivity
Definition: FlowProblem.hpp:1596
ConditionalStorage< enableDispersion, Scalar > dispersivity
Definition: FlowProblem.hpp:1597
Scalar transmissibility
Definition: FlowProblem.hpp:1598