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