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 if (energyModuleType == EnergyModules::SequentialImplicitThermal) {
246 this->enableDriftCompensationTemp_ = Parameters::Get<Parameters::EnableDriftCompensationTemp>();
247 }
248
249 }
250
251 virtual ~FlowProblem() = default;
252
253 void prefetch(const Element& elem) const
254 { this->pffDofData_.prefetch(elem); }
255
267 template <class Restarter>
268 void deserialize(Restarter& res)
269 {
270 // reload the current episode/report step from the deck
271 this->beginEpisode();
272
273 // deserialize the wells
274 wellModel_.deserialize(res);
275
276 // deserialize the aquifer
277 aquiferModel_.deserialize(res);
278 }
279
286 template <class Restarter>
287 void serialize(Restarter& res)
288 {
289 wellModel_.serialize(res);
290
291 aquiferModel_.serialize(res);
292 }
293
294 int episodeIndex() const
295 {
296 return std::max(this->simulator().episodeIndex(), 0);
297 }
298
302 virtual void beginEpisode()
303 {
304 OPM_TIMEBLOCK(beginEpisode);
305 // Proceed to the next report step
306 auto& simulator = this->simulator();
307 int episodeIdx = simulator.episodeIndex();
308 auto& eclState = simulator.vanguard().eclState();
309 const auto& schedule = simulator.vanguard().schedule();
310 const auto& events = schedule[episodeIdx].events();
311
312 if (episodeIdx >= 0 && events.hasEvent(ScheduleEvents::GEO_MODIFIER)) {
313 // bring the contents of the keywords to the current state of the SCHEDULE
314 // section.
315 //
316 // TODO (?): make grid topology changes possible (depending on what exactly
317 // has changed, the grid may need be re-created which has some serious
318 // implications on e.g., the solution of the simulation.)
319 const auto& miniDeck = schedule[episodeIdx].geo_keywords();
320 const auto& cc = simulator.vanguard().grid().comm();
321 eclState.apply_schedule_keywords( miniDeck );
322 eclBroadcast(cc, eclState.getTransMult() );
323
324 // Re-ordering in case of ALUGrid
325 std::function<unsigned int(unsigned int)> equilGridToGrid = [&simulator](unsigned int i) {
326 return simulator.vanguard().gridEquilIdxToGridIdx(i);
327 };
328
329 // re-compute all quantities which may possibly be affected.
330 using TransUpdateQuantities = typename Vanguard::TransmissibilityType::TransUpdateQuantities;
331 transmissibilities_.update(true, TransUpdateQuantities::All, equilGridToGrid);
332 this->referencePorosity_[1] = this->referencePorosity_[0];
334 this->rockFraction_[1] = this->rockFraction_[0];
337 this->model().linearizer().updateDiscretizationParameters();
338 }
339
340 bool tuningEvent = this->beginEpisode_(enableExperiments, this->episodeIndex());
341
342 // set up the wells for the next episode.
343 wellModel_.beginEpisode();
344
345 // set up the aquifers for the next episode.
346 aquiferModel_.beginEpisode();
347
348 // set the size of the initial time step of the episode
349 Scalar dt = limitNextTimeStepSize_(simulator.episodeLength());
350 // negative value of initialTimeStepSize_ indicates no active limit from TSINIT or NEXTSTEP
351 if ( (episodeIdx == 0 || tuningEvent) && this->initialTimeStepSize_ > 0)
352 // allow the size of the initial time step to be set via an external parameter
353 // if TUNING is enabled, also limit the time step size after a tuning event to TSINIT
354 dt = std::min(dt, this->initialTimeStepSize_);
355 simulator.setTimeStepSize(dt);
356 }
357
361 virtual void beginTimeStep()
362 {
363 OPM_TIMEBLOCK(beginTimeStep);
364 const int episodeIdx = this->episodeIndex();
365 const int timeStepSize = this->simulator().timeStepSize();
366
368 episodeIdx,
369 this->simulator().timeStepIndex(),
370 this->simulator().startTime(),
371 this->simulator().time(),
372 timeStepSize,
373 this->simulator().endTime());
374
375 // update maximum water saturation and minimum pressure
376 // used when ROCKCOMP is activated
377 // Do not update max RS first step after a restart
378 this->updateExplicitQuantities_(episodeIdx, timeStepSize, first_step_ && (episodeIdx > 0));
379 first_step_ = false;
380
382 this->model().linearizer().updateBoundaryConditionData();
383 }
384
385 wellModel_.beginTimeStep();
386 aquiferModel_.beginTimeStep();
387 tracerModel_.beginTimeStep();
388 temperatureModel_.beginTimeStep();
389
390 }
391
396 {
397 OPM_TIMEBLOCK(beginIteration);
398 wellModel_.beginIteration();
399 aquiferModel_.beginIteration();
400 }
401
406 {
407 OPM_TIMEBLOCK(endIteration);
408 wellModel_.endIteration();
409 aquiferModel_.endIteration();
410 }
411
415 virtual void endTimeStep()
416 {
417 OPM_TIMEBLOCK(endTimeStep);
418
419#ifndef NDEBUG
420 if constexpr (getPropValue<TypeTag, Properties::EnableDebuggingChecks>()) {
421 // in debug mode, we don't care about performance, so we check
422 // if the model does the right thing (i.e., the mass change
423 // inside the whole reservoir must be equivalent to the fluxes
424 // over the grid's boundaries plus the source rates specified by
425 // the problem).
426 const int rank = this->simulator().gridView().comm().rank();
427 if (rank == 0) {
428 std::cout << "checking conservativeness of solution\n";
429 }
430
431 this->model().checkConservativeness(/*tolerance=*/-1, /*verbose=*/true);
432 if (rank == 0) {
433 std::cout << "solution is sufficiently conservative\n";
434 }
435 }
436#endif // NDEBUG
437
438 auto& simulator = this->simulator();
439 simulator.setTimeStepIndex(simulator.timeStepIndex()+1);
440
441 this->wellModel_.endTimeStep();
442 this->aquiferModel_.endTimeStep();
443 this->tracerModel_.endTimeStep();
444
445 // Compute flux for output
446 this->model().linearizer().updateFlowsInfo();
447
449 OPM_TIMEBLOCK(driftCompansation);
450
451 const auto& residual = this->model().linearizer().residual();
452
453 for (unsigned globalDofIdx = 0; globalDofIdx < residual.size(); globalDofIdx ++) {
454 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(globalDofIdx);
455 this->drift_[sfcdofIdx] = residual[sfcdofIdx] * simulator.timeStepSize();
456
457 if constexpr (getPropValue<TypeTag, Properties::UseVolumetricResidual>()) {
458 this->drift_[sfcdofIdx] *= this->model().dofTotalVolume(sfcdofIdx);
459 }
460 }
461 }
462
463 // Drift compensation needs to be updated before calling the temperature equation
464 if constexpr(energyModuleType == EnergyModules::SequentialImplicitThermal) {
465 this->temperatureModel_.endTimeStep(wellModel_.wellState());
466 }
467 }
468
472 virtual void endEpisode()
473 {
474 const int episodeIdx = this->episodeIndex();
475
476 this->wellModel_.endEpisode();
477 this->aquiferModel_.endEpisode();
478
479 const auto& schedule = this->simulator().vanguard().schedule();
480
481 // End simulation when completed.
482 if (episodeIdx + 1 >= static_cast<int>(schedule.size()) - 1) {
483 this->simulator().setFinished(true);
484 return;
485 }
486
487 // Otherwise, start next episode (report step).
488 this->simulator().startNextEpisode(schedule.stepLength(episodeIdx + 1));
489 }
490
495 virtual void writeOutput(bool verbose)
496 {
497 OPM_TIMEBLOCK(problemWriteOutput);
498
499 if (Parameters::Get<Parameters::EnableWriteAllSolutions>() ||
500 this->episodeWillBeOver())
501 {
502 // Create VTK output as needed.
503 ParentType::writeOutput(verbose);
504 }
505 }
506
510 template <class Context>
511 const DimMatrix& intrinsicPermeability(const Context& context,
512 unsigned spaceIdx,
513 unsigned timeIdx) const
514 {
515 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
516 return transmissibilities_.permeability(globalSpaceIdx);
517 }
518
525 const DimMatrix& intrinsicPermeability(unsigned globalElemIdx) const
526 { return transmissibilities_.permeability(globalElemIdx); }
527
531 template <class Context>
532 Scalar transmissibility(const Context& context,
533 [[maybe_unused]] unsigned fromDofLocalIdx,
534 unsigned toDofLocalIdx) const
535 {
536 assert(fromDofLocalIdx == 0);
537 return pffDofData_.get(context.element(), toDofLocalIdx).transmissibility;
538 }
539
543 Scalar transmissibility(unsigned globalCenterElemIdx, unsigned globalElemIdx) const
544 {
545 return transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
546 }
547
551 template <class Context>
552 Scalar diffusivity(const Context& context,
553 [[maybe_unused]] unsigned fromDofLocalIdx,
554 unsigned toDofLocalIdx) const
555 {
556 assert(fromDofLocalIdx == 0);
557 return *pffDofData_.get(context.element(), toDofLocalIdx).diffusivity;
558 }
559
563 Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
564 return transmissibilities_.diffusivity(globalCellIn, globalCellOut);
565 }
566
570 Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const{
571 return transmissibilities_.dispersivity(globalCellIn, globalCellOut);
572 }
573
577 Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx,
578 const unsigned boundaryFaceIdx) const
579 {
580 return transmissibilities_.thermalTransmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
581 }
582
583
584
585
589 template <class Context>
590 Scalar transmissibilityBoundary(const Context& elemCtx,
591 unsigned boundaryFaceIdx) const
592 {
593 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
594 return transmissibilities_.transmissibilityBoundary(elemIdx, boundaryFaceIdx);
595 }
596
600 Scalar transmissibilityBoundary(const unsigned globalSpaceIdx,
601 const unsigned boundaryFaceIdx) const
602 {
603 return transmissibilities_.transmissibilityBoundary(globalSpaceIdx, boundaryFaceIdx);
604 }
605
606
610 Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn,
611 const unsigned globalSpaceIdxOut) const
612 {
613 return transmissibilities_.thermalHalfTrans(globalSpaceIdxIn,globalSpaceIdxOut);
614 }
615
619 template <class Context>
620 Scalar thermalHalfTransmissibilityIn(const Context& context,
621 unsigned faceIdx,
622 unsigned timeIdx) const
623 {
624 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
625 unsigned toDofLocalIdx = face.exteriorIndex();
626 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransIn;
627 }
628
632 template <class Context>
634 unsigned faceIdx,
635 unsigned timeIdx) const
636 {
637 const auto& face = context.stencil(timeIdx).interiorFace(faceIdx);
638 unsigned toDofLocalIdx = face.exteriorIndex();
639 return *pffDofData_.get(context.element(), toDofLocalIdx).thermalHalfTransOut;
640 }
641
645 template <class Context>
647 unsigned boundaryFaceIdx) const
648 {
649 unsigned elemIdx = elemCtx.globalSpaceIndex(/*dofIdx=*/0, /*timeIdx=*/0);
650 return transmissibilities_.thermalHalfTransBoundary(elemIdx, boundaryFaceIdx);
651 }
652
656 const typename Vanguard::TransmissibilityType& eclTransmissibilities() const
657 { return transmissibilities_; }
658
659
661 { return tracerModel_; }
662
664 { return tracerModel_; }
665
666 TemperatureModel& temperatureModel() // need for restart
667 { return temperatureModel_; }
668
677 template <class Context>
678 Scalar porosity(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
679 {
680 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
681 return this->porosity(globalSpaceIdx, timeIdx);
682 }
683
690 template <class Context>
691 Scalar dofCenterDepth(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
692 {
693 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
694 return this->dofCenterDepth(globalSpaceIdx);
695 }
696
703 Scalar dofCenterDepth(unsigned globalSpaceIdx) const
704 {
705 return this->simulator().vanguard().cellCenterDepth(globalSpaceIdx);
706 }
707
711 template <class Context>
712 Scalar rockCompressibility(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
713 {
714 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
715 return this->rockCompressibility(globalSpaceIdx);
716 }
717
721 template <class Context>
722 Scalar rockReferencePressure(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
723 {
724 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
725 return rockReferencePressure(globalSpaceIdx);
726 }
727
731 Scalar rockReferencePressure(unsigned globalSpaceIdx) const
732 {
733 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
734 if (rock_config.store()) {
735 return asImp_().initialFluidState(globalSpaceIdx).pressure(refPressurePhaseIdx_());
736 }
737 else {
738 if (this->rockParams_.empty())
739 return 1e5;
740
741 unsigned tableIdx = 0;
742 if (!this->rockTableIdx_.empty()) {
743 tableIdx = this->rockTableIdx_[globalSpaceIdx];
744 }
745 return this->rockParams_[tableIdx].referencePressure;
746 }
747 }
748
752 template <class Context>
753 const MaterialLawParams& materialLawParams(const Context& context,
754 unsigned spaceIdx, unsigned timeIdx) const
755 {
756 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
757 return this->materialLawParams(globalSpaceIdx);
758 }
759
760 const MaterialLawParams& materialLawParams(unsigned globalDofIdx) const
761 {
762 return materialLawManager_->materialLawParams(globalDofIdx);
763 }
764
765 const MaterialLawParams& materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
766 {
767 return materialLawManager_->materialLawParams(globalDofIdx, facedir);
768 }
769
773 template <class Context>
775 solidEnergyLawParams(const Context& context,
776 unsigned spaceIdx,
777 unsigned timeIdx) const
778 {
779 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
780 return thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
781 }
782
786 template <class Context>
788 thermalConductionLawParams(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
789 {
790 unsigned globalSpaceIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
791 return thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
792 }
793
800 std::shared_ptr<const EclMaterialLawManager> materialLawManager() const
801 { return materialLawManager_; }
802
803 template <class FluidState, class ...Args>
805 std::array<Evaluation,numPhases> &mobility,
807 FluidState &fluidState,
808 unsigned globalSpaceIdx) const
809 {
810 using ContainerT = std::array<Evaluation, numPhases>;
811 OPM_TIMEBLOCK_LOCAL(updateRelperms, Subsystem::SatProps);
812 {
813 // calculate relative permeabilities. note that we store the result into the
814 // mobility_ class attribute. the division by the phase viscosity happens later.
815 const auto& materialParams = materialLawParams(globalSpaceIdx);
816 MaterialLaw::template relativePermeabilities<ContainerT, FluidState, Args...>(mobility, materialParams, fluidState);
817 Valgrind::CheckDefined(mobility);
818 }
819 if (materialLawManager_->hasDirectionalRelperms()
820 || materialLawManager_->hasDirectionalImbnum())
821 {
822 using Dir = FaceDir::DirEnum;
823 constexpr int ndim = 3;
824 dirMob = std::make_unique<DirectionalMobility<TypeTag>>();
825 Dir facedirs[ndim] = {Dir::XPlus, Dir::YPlus, Dir::ZPlus};
826 for (int i = 0; i<ndim; i++) {
827 const auto& materialParams = materialLawParams(globalSpaceIdx, facedirs[i]);
828 auto& mob_array = dirMob->getArray(i);
829 MaterialLaw::template relativePermeabilities<ContainerT, FluidState, Args...>(mob_array, materialParams, fluidState);
830 }
831 }
832 }
833
837 std::shared_ptr<EclMaterialLawManager> materialLawManager()
838 { return materialLawManager_; }
839
844 template <class Context>
845 unsigned pvtRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
846 { return pvtRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
847
852 template <class Context>
853 unsigned satnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
854 { return this->satnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
855
860 template <class Context>
861 unsigned miscnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
862 { return this->miscnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
863
868 template <class Context>
869 unsigned plmixnumRegionIndex(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
870 { return this->plmixnumRegionIndex(context.globalSpaceIndex(spaceIdx, timeIdx)); }
871
872 // TODO: polymer related might need to go to the blackoil side
877 template <class Context>
878 Scalar maxPolymerAdsorption(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
879 { return this->maxPolymerAdsorption(context.globalSpaceIndex(spaceIdx, timeIdx)); }
880
884 std::string name() const
885 { return this->simulator().vanguard().caseName(); }
886
890 template <class Context>
891 Scalar temperature(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
892 {
893 // use the initial temperature of the DOF if temperature is not a primary
894 // variable
895 unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
896 if constexpr (energyModuleType == EnergyModules::SequentialImplicitThermal)
897 return temperatureModel_.temperature(globalDofIdx);
898
899 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
900 }
901
902
903 Scalar temperature(unsigned globalDofIdx, unsigned /*timeIdx*/) const
904 {
905 // use the initial temperature of the DOF if temperature is not a primary
906 // variable
907 if constexpr (energyModuleType == EnergyModules::SequentialImplicitThermal)
908 return temperatureModel_.temperature(globalDofIdx);
909
910 return asImp_().initialFluidState(globalDofIdx).temperature(/*phaseIdx=*/0);
911 }
912
914 solidEnergyLawParams(unsigned globalSpaceIdx,
915 unsigned /*timeIdx*/) const
916 {
917 return this->thermalLawManager_->solidEnergyLawParams(globalSpaceIdx);
918 }
920 thermalConductionLawParams(unsigned globalSpaceIdx,
921 unsigned /*timeIdx*/)const
922 {
923 return this->thermalLawManager_->thermalConductionLawParams(globalSpaceIdx);
924 }
925
935 Scalar maxOilSaturation(unsigned globalDofIdx) const
936 {
937 if (!this->vapparsActive(this->episodeIndex()))
938 return 0.0;
939
940 return this->maxOilSaturation_[globalDofIdx];
941 }
942
952 void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
953 {
954 if (!this->vapparsActive(this->episodeIndex()))
955 return;
956
957 this->maxOilSaturation_[globalDofIdx] = value;
958 }
959
964 {
965 // Calculate all intensive quantities.
966 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/0);
967
968 // We also need the intensive quantities for timeIdx == 1
969 // corresponding to the start of the current timestep, if we
970 // do not use the storage cache, or if we cannot recycle the
971 // first iteration storage.
972 if (!this->model().enableStorageCache() || !this->recycleFirstIterationStorage()) {
973 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx*/1);
974 }
975
976 // initialize the wells. Note that this needs to be done after initializing the
977 // intrinsic permeabilities and the after applying the initial solution because
978 // the well model uses these...
979 wellModel_.init();
980
981 aquiferModel_.initialSolutionApplied();
982
983 const bool invalidateFromHyst = updateHysteresis_();
984 if (invalidateFromHyst) {
985 OPM_TIMEBLOCK(beginTimeStepInvalidateIntensiveQuantities);
986 this->model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0);
987 }
988 }
989
995 template <class Context>
996 void source(RateVector& rate,
997 const Context& context,
998 unsigned spaceIdx,
999 unsigned timeIdx) const
1000 {
1001 const unsigned globalDofIdx = context.globalSpaceIndex(spaceIdx, timeIdx);
1002 source(rate, globalDofIdx, timeIdx);
1003 }
1004
1005 void source(RateVector& rate,
1006 unsigned globalDofIdx,
1007 unsigned timeIdx) const
1008 {
1009 OPM_TIMEBLOCK_LOCAL(eclProblemSource, Subsystem::Assembly);
1010 rate = 0.0;
1011
1012 // Add well contribution to source here.
1013 wellModel_.computeTotalRatesForDof(rate, globalDofIdx);
1014
1015 // convert the source term from the total mass rate of the
1016 // cell to the one per unit of volume as used by the model.
1017 for (unsigned eqIdx = 0; eqIdx < numEq; ++ eqIdx) {
1018 rate[eqIdx] /= this->model().dofTotalVolume(globalDofIdx);
1019
1020 Valgrind::CheckDefined(rate[eqIdx]);
1021 assert(isfinite(rate[eqIdx]));
1022 }
1023
1024 // Add non-well sources.
1025 addToSourceDense(rate, globalDofIdx, timeIdx);
1026 }
1027
1028 virtual void addToSourceDense(RateVector& rate,
1029 unsigned globalDofIdx,
1030 unsigned timeIdx) const = 0;
1031
1037 const WellModel& wellModel() const
1038 { return wellModel_; }
1039
1041 { return wellModel_; }
1042
1044 { return aquiferModel_; }
1045
1047 { return aquiferModel_; }
1048
1051
1059 {
1060 OPM_TIMEBLOCK(nexTimeStepSize);
1061 // allow external code to do the timestepping
1062 if (this->nextTimeStepSize_ > 0.0)
1063 return this->nextTimeStepSize_;
1064
1065 const auto& simulator = this->simulator();
1066 int episodeIdx = simulator.episodeIndex();
1067
1068 // for the initial episode, we use a fixed time step size
1069 if (episodeIdx < 0)
1070 return this->initialTimeStepSize_;
1071
1072 // ask the newton method for a suggestion. This suggestion will be based on how
1073 // well the previous time step converged. After that, apply the runtime time
1074 // stepping constraints.
1075 const auto& newtonMethod = this->model().newtonMethod();
1076 return limitNextTimeStepSize_(newtonMethod.suggestTimeStepSize(simulator.timeStepSize()));
1077 }
1078
1084 template <class LhsEval>
1085 LhsEval rockCompPoroMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1086 {
1087 OPM_TIMEBLOCK_LOCAL(rockCompPoroMultiplier, Subsystem::PvtProps);
1088 if (this->rockCompPoroMult_.empty() && this->rockCompPoroMultWc_.empty())
1089 return 1.0;
1090
1091 unsigned tableIdx = 0;
1092 if (!this->rockTableIdx_.empty())
1093 tableIdx = this->rockTableIdx_[elementIdx];
1094
1095 const auto& fs = intQuants.fluidState();
1096 LhsEval effectivePressure = decay<LhsEval>(fs.pressure(refPressurePhaseIdx_()));
1097 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1098 if (!this->minRefPressure_.empty())
1099 // The pore space change is irreversible
1100 effectivePressure =
1101 min(decay<LhsEval>(fs.pressure(refPressurePhaseIdx_())),
1102 this->minRefPressure_[elementIdx]);
1103
1104 if (!this->overburdenPressure_.empty())
1105 effectivePressure -= this->overburdenPressure_[elementIdx];
1106
1107 if (rock_config.store()) {
1108 effectivePressure -= asImp_().initialFluidState(elementIdx).pressure(refPressurePhaseIdx_());
1109 }
1110
1111 if (!this->rockCompPoroMult_.empty()) {
1112 return this->rockCompPoroMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1113 }
1114
1115 // water compaction
1116 assert(!this->rockCompPoroMultWc_.empty());
1117 LhsEval SwMax = max(decay<LhsEval>(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1118 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1119
1120 return this->rockCompPoroMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1121 }
1122
1128 template <class LhsEval>
1129 LhsEval rockCompTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1130 {
1131 auto obtain = [](const auto& value)
1132 {
1133 if constexpr (std::is_same_v<LhsEval, Scalar>) {
1134 return getValue(value);
1135 } else {
1136 return value;
1137 }
1138 };
1139 return rockCompTransMultiplier<LhsEval>(intQuants, elementIdx, obtain);
1140 }
1141
1142 template <class LhsEval, class Callback>
1143 LhsEval rockCompTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx, Callback& obtain) const
1144 {
1145 const bool implicit = !this->explicitRockCompaction_;
1146 return implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx, obtain)
1147 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1148 }
1149
1150 template <class LhsEval, class Callback>
1151 LhsEval wellTransMultiplier(const IntensiveQuantities& intQuants, unsigned elementIdx, Callback& obtain) const
1152 {
1153 OPM_TIMEBLOCK_LOCAL(wellTransMultiplier, Subsystem::Wells);
1154
1155 const bool implicit = !this->explicitRockCompaction_;
1156 LhsEval trans_mult = implicit ? this->simulator().problem().template computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx, obtain)
1157 : this->simulator().problem().getRockCompTransMultVal(elementIdx);
1158 trans_mult *= this->simulator().problem().template permFactTransMultiplier<LhsEval>(intQuants, elementIdx, obtain);
1159
1160 return trans_mult;
1161 }
1162
1163 std::pair<BCType, RateVector> boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
1164 {
1165 OPM_TIMEBLOCK_LOCAL(boundaryCondition, Subsystem::Assembly);
1167 return { BCType::NONE, RateVector(0.0) };
1168 }
1169 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
1170 const auto& schedule = this->simulator().vanguard().schedule();
1171 if (bcindex_(dir)[globalSpaceIdx] == 0) {
1172 return { BCType::NONE, RateVector(0.0) };
1173 }
1174 if (schedule[this->episodeIndex()].bcprop.size() == 0) {
1175 return { BCType::NONE, RateVector(0.0) };
1176 }
1177 const auto& bc = schedule[this->episodeIndex()].bcprop[bcindex_(dir)[globalSpaceIdx]];
1178 if (bc.bctype!=BCType::RATE) {
1179 return { bc.bctype, RateVector(0.0) };
1180 }
1181
1182 RateVector rate = 0.0;
1183 switch (bc.component) {
1184 case BCComponent::OIL:
1185 rate[FluidSystem::canonicalToActiveCompIdx(oilCompIdx)] = bc.rate;
1186 break;
1187 case BCComponent::GAS:
1188 rate[FluidSystem::canonicalToActiveCompIdx(gasCompIdx)] = bc.rate;
1189 break;
1190 case BCComponent::WATER:
1191 rate[FluidSystem::canonicalToActiveCompIdx(waterCompIdx)] = bc.rate;
1192 break;
1193 case BCComponent::SOLVENT:
1194 this->handleSolventBC(bc, rate);
1195 break;
1196 case BCComponent::POLYMER:
1197 this->handlePolymerBC(bc, rate);
1198 break;
1199 case BCComponent::MICR:
1200 this->handleMicrBC(bc, rate);
1201 break;
1202 case BCComponent::OXYG:
1203 this->handleOxygBC(bc, rate);
1204 break;
1205 case BCComponent::UREA:
1206 this->handleUreaBC(bc, rate);
1207 break;
1208 case BCComponent::NONE:
1209 throw std::logic_error("you need to specify the component when RATE type is set in BC");
1210 break;
1211 }
1212 //TODO add support for enthalpy rate
1213 return {bc.bctype, rate};
1214 }
1215
1216
1217 template<class Serializer>
1218 void serializeOp(Serializer& serializer)
1219 {
1220 serializer(static_cast<BaseType&>(*this));
1221 serializer(drift_);
1222 serializer(wellModel_);
1223 serializer(aquiferModel_);
1224 serializer(tracerModel_);
1225 serializer(*materialLawManager_);
1226 }
1227
1228 const GlobalEqVector& drift() const
1229 {
1230 return drift_;
1231 }
1232
1233private:
1234 Implementation& asImp_()
1235 { return *static_cast<Implementation *>(this); }
1236
1237 const Implementation& asImp_() const
1238 { return *static_cast<const Implementation *>(this); }
1239
1240protected:
1241 template<class UpdateFunc>
1242 void updateProperty_(const std::string& failureMsg,
1243 UpdateFunc func)
1244 {
1245 OPM_TIMEBLOCK(updateProperty);
1246 const auto& model = this->simulator().model();
1247 const auto& primaryVars = model.solution(/*timeIdx*/0);
1248 const auto& vanguard = this->simulator().vanguard();
1249 std::size_t numGridDof = primaryVars.size();
1251#ifdef _OPENMP
1252#pragma omp parallel for
1253#endif
1254 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++dofIdx) {
1255 const auto& iq = *model.cachedIntensiveQuantities(dofIdx, /*timeIdx=*/ 0);
1256 func(dofIdx, iq);
1257 }
1258 OPM_END_PARALLEL_TRY_CATCH(failureMsg, vanguard.grid().comm());
1259 }
1260
1262 {
1263 OPM_TIMEBLOCK(updateMaxOilSaturation);
1264 int episodeIdx = this->episodeIndex();
1265
1266 // we use VAPPARS
1267 if (this->vapparsActive(episodeIdx)) {
1268 this->updateProperty_("FlowProblem::updateMaxOilSaturation_() failed:",
1269 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1270 {
1271 this->updateMaxOilSaturation_(compressedDofIdx,iq);
1272 });
1273 return true;
1274 }
1275
1276 return false;
1277 }
1278
1279 bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1280 {
1281 OPM_TIMEBLOCK_LOCAL(updateMaxOilSaturation, Subsystem::SatProps);
1282 const auto& fs = iq.fluidState();
1283 const Scalar So = decay<Scalar>(fs.saturation(refPressurePhaseIdx_()));
1284 auto& mos = this->maxOilSaturation_;
1285 if(mos[compressedDofIdx] < So){
1286 mos[compressedDofIdx] = So;
1287 return true;
1288 }else{
1289 return false;
1290 }
1291 }
1292
1294 {
1295 OPM_TIMEBLOCK(updateMaxWaterSaturation);
1296 // water compaction is activated in ROCKCOMP
1297 if (this->maxWaterSaturation_.empty())
1298 return false;
1299
1300 this->maxWaterSaturation_[/*timeIdx=*/1] = this->maxWaterSaturation_[/*timeIdx=*/0];
1301 this->updateProperty_("FlowProblem::updateMaxWaterSaturation_() failed:",
1302 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1303 {
1304 this->updateMaxWaterSaturation_(compressedDofIdx,iq);
1305 });
1306 return true;
1307 }
1308
1309
1310 bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1311 {
1312 OPM_TIMEBLOCK_LOCAL(updateMaxWaterSaturation, Subsystem::SatProps);
1313 const auto& fs = iq.fluidState();
1314 const Scalar Sw = decay<Scalar>(fs.saturation(waterPhaseIdx));
1315 auto& mow = this->maxWaterSaturation_;
1316 if(mow[compressedDofIdx]< Sw){
1317 mow[compressedDofIdx] = Sw;
1318 return true;
1319 }else{
1320 return false;
1321 }
1322 }
1323
1325 {
1326 OPM_TIMEBLOCK(updateMinPressure);
1327 // IRREVERS option is used in ROCKCOMP
1328 if (this->minRefPressure_.empty())
1329 return false;
1330
1331 this->updateProperty_("FlowProblem::updateMinPressure_() failed:",
1332 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1333 {
1334 this->updateMinPressure_(compressedDofIdx,iq);
1335 });
1336 return true;
1337 }
1338
1339 bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities& iq){
1340 OPM_TIMEBLOCK_LOCAL(updateMinPressure, Subsystem::PvtProps);
1341 const auto& fs = iq.fluidState();
1342 const Scalar min_pressure = getValue(fs.pressure(refPressurePhaseIdx_()));
1343 auto& min_pressures = this->minRefPressure_;
1344 if(min_pressures[compressedDofIdx]> min_pressure){
1345 min_pressures[compressedDofIdx] = min_pressure;
1346 return true;
1347 }else{
1348 return false;
1349 }
1350 }
1351
1352 // \brief Function to assign field properties of type double, 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 std::function<std::vector<double>(const FieldPropsManager&, const std::string&)>
1358 {
1359 const auto& lookup = this->lookUpData_;
1360 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString)
1361 {
1362 return lookup.assignFieldPropsDoubleOnLeaf(fieldPropManager, propString);
1363 };
1364 }
1365
1366 // \brief Function to assign field properties of type int, unsigned int, ..., on the leaf grid view.
1367 //
1368 // For CpGrid with local grid refinement, the field property of a cell on the leaf
1369 // is inherited from its parent or equivalent (when has no parent) cell on level zero.
1370 template<typename IntType>
1371 std::function<std::vector<IntType>(const FieldPropsManager&, const std::string&, bool)>
1373 {
1374 const auto& lookup = this->lookUpData_;
1375 return [&lookup](const FieldPropsManager& fieldPropManager, const std::string& propString, bool needsTranslation)
1376 {
1377 return lookup.template assignFieldPropsIntOnLeaf<IntType>(fieldPropManager, propString, needsTranslation);
1378 };
1379 }
1380
1382 {
1383 OPM_TIMEBLOCK(readMaterialParameters);
1384 const auto& simulator = this->simulator();
1385 const auto& vanguard = simulator.vanguard();
1386 const auto& eclState = vanguard.eclState();
1387
1388 // the PVT and saturation region numbers
1390 this->updatePvtnum_();
1391 this->updateSatnum_();
1392
1393 // the MISC region numbers (solvent model)
1394 this->updateMiscnum_();
1395 // the PLMIX region numbers (polymer model)
1396 this->updatePlmixnum_();
1397
1398 OPM_END_PARALLEL_TRY_CATCH("Invalid region numbers: ", vanguard.gridView().comm());
1400 // porosity
1402 this->referencePorosity_[1] = this->referencePorosity_[0];
1404
1406 // rock fraction
1408 this->rockFraction_[1] = this->rockFraction_[0];
1410
1412 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1413 materialLawManager_ = std::make_shared<EclMaterialLawManager>();
1414 materialLawManager_->initFromState(eclState);
1415 materialLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1416 this-> template fieldPropIntTypeOnLeafAssigner_<int>(),
1419 }
1420
1422 {
1423 if constexpr (energyModuleType == EnergyModules::FullyImplicitThermal ||
1424 energyModuleType == EnergyModules::SequentialImplicitThermal )
1425 {
1426 const auto& simulator = this->simulator();
1427 const auto& vanguard = simulator.vanguard();
1428 const auto& eclState = vanguard.eclState();
1429
1430 // fluid-matrix interactions (saturation functions; relperm/capillary pressure)
1431 thermalLawManager_ = std::make_shared<EclThermalLawManager>();
1432 thermalLawManager_->initParamsForElements(eclState, this->model().numGridDof(),
1434 this-> template fieldPropIntTypeOnLeafAssigner_<unsigned int>());
1435 }
1436 }
1437
1439 {
1440 const auto& simulator = this->simulator();
1441 const auto& vanguard = simulator.vanguard();
1442 const auto& eclState = vanguard.eclState();
1443
1444 std::size_t numDof = this->model().numGridDof();
1445
1446 this->referencePorosity_[/*timeIdx=*/0].resize(numDof);
1447
1448 const auto& fp = eclState.fieldProps();
1449 const std::vector<double> porvData = this -> fieldPropDoubleOnLeafAssigner_()(fp, "PORV");
1450 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1451 int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(dofIdx);
1452 Scalar poreVolume = porvData[dofIdx];
1453
1454 // we define the porosity as the accumulated pore volume divided by the
1455 // geometric volume of the element. Note that -- in pathetic cases -- it can
1456 // be larger than 1.0!
1457 Scalar dofVolume = simulator.model().dofTotalVolume(sfcdofIdx);
1458 assert(dofVolume > 0.0);
1459 this->referencePorosity_[/*timeIdx=*/0][sfcdofIdx] = poreVolume/dofVolume;
1460 }
1461 }
1462
1464
1465 const bool solveEnergyEquation = (energyModuleType == EnergyModules::FullyImplicitThermal ||
1466 energyModuleType == EnergyModules::SequentialImplicitThermal);
1467 if (!solveEnergyEquation)
1468 return;
1469
1470 const auto& simulator = this->simulator();
1471 const auto& vanguard = simulator.vanguard();
1472 const auto& eclState = vanguard.eclState();
1473
1474 std::size_t numDof = this->model().numGridDof();
1475 this->rockFraction_[/*timeIdx=*/0].resize(numDof);
1476 // For the energy equation, we need the volume of the rock.
1477 // The volume of the rock is computed by rockFraction * geometric volume of the element.
1478 // The reference porosity is defined as porosity * ntg * pore-volume-multiplier.
1479 // A common practice in reservoir simulation is to use large pore-volume-multipliers in boundary cells
1480 // to model boundary conditions other than no-flow. This may result in reference porosities that are larger than 1.
1481 // A simple (1-reference porosity) * geometric volume of the element may give unphysical results.
1482 // We therefore instead consider the pore-volume-multiplier as a volume multiplier. The rock fraction is thus given by
1483 // (1 - porosity * ntg) * pore-volume-multiplier = (1 - porosity * ntg) * reference porosity / (porosity * ntg)
1484 const auto& fp = eclState.fieldProps();
1485 const std::vector<double> poroData = this->fieldPropDoubleOnLeafAssigner_()(fp, "PORO");
1486 const std::vector<double> ntgData = this->fieldPropDoubleOnLeafAssigner_()(fp, "NTG");
1487
1488 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
1489 const auto ntg = ntgData[dofIdx];
1490 const auto poro_eff = ntg * poroData[dofIdx];
1491 const int sfcdofIdx = simulator.vanguard().gridEquilIdxToGridIdx(dofIdx);
1492 const auto rock_fraction = (1 - poro_eff) * this->referencePorosity_[/*timeIdx=*/0][sfcdofIdx] / poro_eff;
1493 this->rockFraction_[/*timeIdx=*/0][sfcdofIdx] = rock_fraction;
1494 }
1495 }
1496
1498 {
1499 // TODO: whether we should move this to FlowProblemBlackoil
1500 const auto& simulator = this->simulator();
1501 const auto& vanguard = simulator.vanguard();
1502 const auto& eclState = vanguard.eclState();
1503
1504 if (eclState.getInitConfig().hasEquil())
1506 else
1508
1509 //initialize min/max values
1510 std::size_t numElems = this->model().numGridDof();
1511 for (std::size_t elemIdx = 0; elemIdx < numElems; ++elemIdx) {
1512 const auto& fs = asImp_().initialFluidStates()[elemIdx];
1513 if (!this->maxWaterSaturation_.empty() && waterPhaseIdx > -1)
1514 this->maxWaterSaturation_[elemIdx] = std::max(this->maxWaterSaturation_[elemIdx], fs.saturation(waterPhaseIdx));
1515 if (!this->maxOilSaturation_.empty() && oilPhaseIdx > -1)
1516 this->maxOilSaturation_[elemIdx] = std::max(this->maxOilSaturation_[elemIdx], fs.saturation(oilPhaseIdx));
1517 if (!this->minRefPressure_.empty() && refPressurePhaseIdx_() > -1)
1518 this->minRefPressure_[elemIdx] = std::min(this->minRefPressure_[elemIdx], fs.pressure(refPressurePhaseIdx_()));
1519 }
1520 }
1521
1522 virtual void readEquilInitialCondition_() = 0;
1524
1525 // update the hysteresis parameters of the material laws for the whole grid
1527 {
1528 if (!materialLawManager_->enableHysteresis())
1529 return false;
1530
1531 // we need to update the hysteresis data for _all_ elements (i.e., not just the
1532 // interior ones) to avoid desynchronization of the processes in the parallel case!
1533 this->updateProperty_("FlowProblem::updateHysteresis_() failed:",
1534 [this](unsigned compressedDofIdx, const IntensiveQuantities& iq)
1535 {
1536 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1537 });
1538 return true;
1539 }
1540
1541
1542 bool updateHysteresis_(unsigned compressedDofIdx, const IntensiveQuantities& iq)
1543 {
1544 OPM_TIMEBLOCK_LOCAL(updateHysteresis_, Subsystem::SatProps);
1545 materialLawManager_->updateHysteresis(iq.fluidState(), compressedDofIdx);
1546 //TODO change materials to give a bool
1547 return true;
1548 }
1549
1550 Scalar getRockCompTransMultVal(std::size_t dofIdx) const
1551 {
1552 if (this->rockCompTransMultVal_.empty())
1553 return 1.0;
1554
1555 return this->rockCompTransMultVal_[dofIdx];
1556 }
1557
1558protected:
1560 {
1561 ConditionalStorage<enableFullyImplicitThermal, Scalar> thermalHalfTransIn;
1562 ConditionalStorage<enableFullyImplicitThermal, Scalar> thermalHalfTransOut;
1563 ConditionalStorage<enableDiffusion, Scalar> diffusivity;
1564 ConditionalStorage<enableDispersion, Scalar> dispersivity;
1566 };
1567
1568 // update the prefetch friendly data object
1570 {
1571 const auto& distFn =
1572 [this](PffDofData_& dofData,
1573 const Stencil& stencil,
1574 unsigned localDofIdx)
1575 -> void
1576 {
1577 const auto& elementMapper = this->model().elementMapper();
1578
1579 unsigned globalElemIdx = elementMapper.index(stencil.entity(localDofIdx));
1580 if (localDofIdx != 0) {
1581 unsigned globalCenterElemIdx = elementMapper.index(stencil.entity(/*dofIdx=*/0));
1582 dofData.transmissibility = transmissibilities_.transmissibility(globalCenterElemIdx, globalElemIdx);
1583
1584 if constexpr (enableFullyImplicitThermal) {
1585 *dofData.thermalHalfTransIn = transmissibilities_.thermalHalfTrans(globalCenterElemIdx, globalElemIdx);
1586 *dofData.thermalHalfTransOut = transmissibilities_.thermalHalfTrans(globalElemIdx, globalCenterElemIdx);
1587 }
1588 if constexpr (enableDiffusion)
1589 *dofData.diffusivity = transmissibilities_.diffusivity(globalCenterElemIdx, globalElemIdx);
1590 if (enableDispersion)
1591 dofData.dispersivity = transmissibilities_.dispersivity(globalCenterElemIdx, globalElemIdx);
1592 }
1593 };
1594
1595 pffDofData_.update(distFn);
1596 }
1597
1598 virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart) = 0;
1599
1601 {
1602 const auto& simulator = this->simulator();
1603 const auto& vanguard = simulator.vanguard();
1604 const auto& bcconfig = vanguard.eclState().getSimulationConfig().bcconfig();
1605 if (bcconfig.size() > 0) {
1607
1608 std::size_t numCartDof = vanguard.cartesianSize();
1609 unsigned numElems = vanguard.gridView().size(/*codim=*/0);
1610 std::vector<int> cartesianToCompressedElemIdx(numCartDof, -1);
1611
1612 for (unsigned elemIdx = 0; elemIdx < numElems; ++elemIdx)
1613 cartesianToCompressedElemIdx[vanguard.cartesianIndex(elemIdx)] = elemIdx;
1614
1615 bcindex_.resize(numElems, 0);
1616 auto loopAndApply = [&cartesianToCompressedElemIdx,
1617 &vanguard](const auto& bcface,
1618 auto apply)
1619 {
1620 for (int i = bcface.i1; i <= bcface.i2; ++i) {
1621 for (int j = bcface.j1; j <= bcface.j2; ++j) {
1622 for (int k = bcface.k1; k <= bcface.k2; ++k) {
1623 std::array<int, 3> tmp = {i,j,k};
1624 auto elemIdx = cartesianToCompressedElemIdx[vanguard.cartesianIndex(tmp)];
1625 if (elemIdx >= 0)
1626 apply(elemIdx);
1627 }
1628 }
1629 }
1630 };
1631 for (const auto& bcface : bcconfig) {
1632 std::vector<int>& data = bcindex_(bcface.dir);
1633 const int index = bcface.index;
1634 loopAndApply(bcface,
1635 [&data,index](int elemIdx)
1636 { data[elemIdx] = index; });
1637 }
1638 }
1639 }
1640
1641 // this method applies the runtime constraints specified via the deck and/or command
1642 // line parameters for the size of the next time step.
1644 {
1645 if constexpr (enableExperiments) {
1646 const auto& simulator = this->simulator();
1647 const auto& schedule = simulator.vanguard().schedule();
1648 int episodeIdx = simulator.episodeIndex();
1649
1650 // first thing in the morning, limit the time step size to the maximum size
1651 Scalar maxTimeStepSize = Parameters::Get<Parameters::SolverMaxTimeStepInDays<Scalar>>() * 24 * 60 * 60;
1652 int reportStepIdx = std::max(episodeIdx, 0);
1653 if (this->enableTuning_) {
1654 const auto& tuning = schedule[reportStepIdx].tuning();
1655 maxTimeStepSize = tuning.TSMAXZ;
1656 }
1657
1658 dtNext = std::min(dtNext, maxTimeStepSize);
1659
1660 Scalar remainingEpisodeTime =
1661 simulator.episodeStartTime() + simulator.episodeLength()
1662 - (simulator.startTime() + simulator.time());
1663 assert(remainingEpisodeTime >= 0.0);
1664
1665 // if we would have a small amount of time left over in the current episode, make
1666 // two equal time steps instead of a big and a small one
1667 if (remainingEpisodeTime/2.0 < dtNext && dtNext < remainingEpisodeTime*(1.0 - 1e-5))
1668 // note: limiting to the maximum time step size here is probably not strictly
1669 // necessary, but it should not hurt and is more fool-proof
1670 dtNext = std::min(maxTimeStepSize, remainingEpisodeTime/2.0);
1671
1672 if (simulator.episodeStarts()) {
1673 // if a well event occurred, respect the limit for the maximum time step after
1674 // that, too
1675 const auto& events = simulator.vanguard().schedule()[reportStepIdx].events();
1676 bool wellEventOccured =
1677 events.hasEvent(ScheduleEvents::NEW_WELL)
1678 || events.hasEvent(ScheduleEvents::PRODUCTION_UPDATE)
1679 || events.hasEvent(ScheduleEvents::INJECTION_UPDATE)
1680 || events.hasEvent(ScheduleEvents::WELL_STATUS_CHANGE);
1681 if (episodeIdx >= 0 && wellEventOccured && this->maxTimeStepAfterWellEvent_ > 0)
1682 dtNext = std::min(dtNext, this->maxTimeStepAfterWellEvent_);
1683 }
1684 }
1685
1686 return dtNext;
1687 }
1688
1690 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1691 return oilPhaseIdx;
1692 }
1693 else if (FluidSystem::phaseIsActive(gasPhaseIdx)) {
1694 return gasPhaseIdx;
1695 }
1696 else {
1697 return waterPhaseIdx;
1698 }
1699 }
1700
1702 {
1703 const auto& model = this->simulator().model();
1704 std::size_t numGridDof = this->model().numGridDof();
1705 this->rockCompTransMultVal_.resize(numGridDof, 1.0);
1706 for (std::size_t elementIdx = 0; elementIdx < numGridDof; ++elementIdx) {
1707 const auto& iq = *model.cachedIntensiveQuantities(elementIdx, /*timeIdx=*/ 0);
1708 Scalar trans_mult = computeRockCompTransMultiplier_<Scalar>(iq, elementIdx);
1709 this->rockCompTransMultVal_[elementIdx] = trans_mult;
1710 }
1711 }
1712
1718 template <class LhsEval>
1719 LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities& intQuants, unsigned elementIdx) const
1720 {
1721 auto obtain = [](const auto& value)
1722 {
1723 if constexpr (std::is_same_v<LhsEval, Scalar>) {
1724 return getValue(value);
1725 } else {
1726 return value;
1727 }
1728 };
1729
1730 return computeRockCompTransMultiplier_<LhsEval>(intQuants, elementIdx, obtain);
1731 }
1732
1733 template <class LhsEval, class Callback>
1734 LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities& intQuants, unsigned elementIdx, Callback& obtain) const
1735 {
1736 OPM_TIMEBLOCK_LOCAL(computeRockCompTransMultiplier, Subsystem::PvtProps);
1737 if (this->rockCompTransMult_.empty() && this->rockCompTransMultWc_.empty())
1738 return 1.0;
1739
1740 unsigned tableIdx = 0;
1741 if (!this->rockTableIdx_.empty())
1742 tableIdx = this->rockTableIdx_[elementIdx];
1743
1744 const auto& fs = intQuants.fluidState();
1745 LhsEval effectivePressure = obtain(fs.pressure(refPressurePhaseIdx_()));
1746 const auto& rock_config = this->simulator().vanguard().eclState().getSimulationConfig().rock_config();
1747 if (!this->minRefPressure_.empty())
1748 // The pore space change is irreversible
1749 effectivePressure =
1750 min(obtain(fs.pressure(refPressurePhaseIdx_())),
1751 this->minRefPressure_[elementIdx]);
1752
1753 if (!this->overburdenPressure_.empty())
1754 effectivePressure -= this->overburdenPressure_[elementIdx];
1755
1756 if (rock_config.store()) {
1757 effectivePressure -= asImp_().initialFluidState(elementIdx).pressure(refPressurePhaseIdx_());
1758 }
1759
1760 if (!this->rockCompTransMult_.empty())
1761 return this->rockCompTransMult_[tableIdx].eval(effectivePressure, /*extrapolation=*/true);
1762
1763 // water compaction
1764 assert(!this->rockCompTransMultWc_.empty());
1765 LhsEval SwMax = max(obtain(fs.saturation(waterPhaseIdx)), this->maxWaterSaturation_[elementIdx]);
1766 LhsEval SwDeltaMax = SwMax - asImp_().initialFluidStates()[elementIdx].saturation(waterPhaseIdx);
1767
1768 return this->rockCompTransMultWc_[tableIdx].eval(effectivePressure, SwDeltaMax, /*extrapolation=*/true);
1769 }
1770
1771 typename Vanguard::TransmissibilityType transmissibilities_;
1772
1773 std::shared_ptr<EclMaterialLawManager> materialLawManager_;
1774 std::shared_ptr<EclThermalLawManager> thermalLawManager_;
1775
1777
1780
1784
1785 template<class T>
1786 struct BCData
1787 {
1788 std::array<std::vector<T>,6> data;
1789
1790 void resize(std::size_t size, T defVal)
1791 {
1792 for (auto& d : data)
1793 d.resize(size, defVal);
1794 }
1795
1796 const std::vector<T>& operator()(FaceDir::DirEnum dir) const
1797 {
1798 if (dir == FaceDir::DirEnum::Unknown)
1799 throw std::runtime_error("Tried to access BC data for the 'Unknown' direction");
1800 int idx = 0;
1801 int div = static_cast<int>(dir);
1802 while ((div /= 2) >= 1)
1803 ++idx;
1804 assert(idx >= 0 && idx <= 5);
1805 return data[idx];
1806 }
1807
1808 std::vector<T>& operator()(FaceDir::DirEnum dir)
1809 {
1810 return const_cast<std::vector<T>&>(std::as_const(*this)(dir));
1811 }
1812 };
1813
1814 virtual void handleSolventBC(const BCProp::BCFace&, RateVector&) const = 0;
1815
1816 virtual void handlePolymerBC(const BCProp::BCFace&, RateVector&) const = 0;
1817
1818 virtual void handleMicrBC(const BCProp::BCFace&, RateVector&) const = 0;
1819
1820 virtual void handleOxygBC(const BCProp::BCFace&, RateVector&) const = 0;
1821
1822 virtual void handleUreaBC(const BCProp::BCFace&, RateVector&) const = 0;
1823
1826 bool first_step_ = true;
1827
1830 virtual bool episodeWillBeOver() const
1831 {
1832 return this->simulator().episodeWillBeOver();
1833 }
1834};
1835
1836} // namespace Opm
1837
1838#endif // OPM_FLOW_PROBLEM_HPP
#define OPM_END_PARALLEL_TRY_CATCH(prefix, comm)
Catch exception and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:192
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition: FlowGenericProblem.hpp:61
Scalar maxPolymerAdsorption(unsigned elemIdx) const
Returns the max polymer adsorption value.
Definition: FlowGenericProblem_impl.hpp:733
unsigned pvtRegionIndex(unsigned elemIdx) const
Returns the index the relevant PVT region given a cell index.
Definition: FlowGenericProblem_impl.hpp:692
std::array< std::vector< Scalar >, 2 > rockFraction_
Definition: FlowGenericProblem.hpp:323
std::function< unsigned(unsigned)> lookupIdxOnLevelZeroAssigner_()
Definition: FlowGenericProblem.hpp:369
std::vector< TabulatedTwoDFunction > rockCompPoroMultWc_
Definition: FlowGenericProblem.hpp:332
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:712
unsigned satnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant saturation function region given a cell index.
Definition: FlowGenericProblem_impl.hpp:702
void beginTimeStep_(bool enableExperiments, int episodeIdx, int timeStepIndex, Scalar startTime, Scalar time, Scalar timeStepSize, Scalar endTime)
Definition: FlowGenericProblem_impl.hpp:444
std::vector< TabulatedTwoDFunction > rockCompTransMultWc_
Definition: FlowGenericProblem.hpp:333
unsigned plmixnumRegionIndex(unsigned elemIdx) const
Returns the index the relevant PLMIXNUM (for polymer module) region given a cell index.
Definition: FlowGenericProblem_impl.hpp:722
std::array< std::vector< Scalar >, 2 > referencePorosity_
Definition: FlowGenericProblem.hpp:322
bool beginEpisode_(bool enableExperiments, int episodeIdx)
Definition: FlowGenericProblem_impl.hpp:407
bool shouldWriteOutput() const
Always returns true. The ecl output writer takes care of the rest.
Definition: FlowGenericProblem.hpp:279
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:288
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:1830
const WellModel & wellModel() const
Returns a reference to the ECL well manager used by the problem.
Definition: FlowProblem.hpp:1037
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:543
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
@ enableSolvent
Definition: FlowProblem.hpp:134
@ enableBioeffects
Definition: FlowProblem.hpp:120
bool nonTrivialBoundaryConditions() const
Definition: FlowProblem.hpp:1049
const GlobalEqVector & drift() const
Definition: FlowProblem.hpp:1228
virtual void writeOutput(bool verbose)
Write the requested quantities of the current solution into the output files.
Definition: FlowProblem.hpp:495
std::function< std::vector< IntType >(const FieldPropsManager &, const std::string &, bool)> fieldPropIntTypeOnLeafAssigner_()
Definition: FlowProblem.hpp:1372
const DimMatrix & intrinsicPermeability(unsigned globalElemIdx) const
This method returns the intrinsic permeability tensor given a global element index.
Definition: FlowProblem.hpp:525
unsigned pvtRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:845
LhsEval wellTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Definition: FlowProblem.hpp:1151
Scalar porosity(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:678
void beginIteration()
Called by the simulator before each Newton-Raphson iteration.
Definition: FlowProblem.hpp:395
GetPropType< TypeTag, Properties::Vanguard > Vanguard
Definition: FlowProblem.hpp:108
GetPropType< TypeTag, Properties::DofMapper > DofMapper
Definition: FlowProblem.hpp:158
@ oilPhaseIdx
Definition: FlowProblem.hpp:138
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:1542
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:853
virtual void handleOxygBC(const BCProp::BCFace &, RateVector &) const =0
virtual void updateExplicitQuantities_(int episodeIdx, int timeStepSize, bool first_step_after_restart)=0
@ enablePolymerMolarWeight
Definition: FlowProblem.hpp:132
virtual void handleUreaBC(const BCProp::BCFace &, RateVector &) const =0
Scalar thermalHalfTransmissibilityOut(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:633
@ gasCompIdx
Definition: FlowProblem.hpp:143
bool first_step_
Definition: FlowProblem.hpp:1826
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: FlowProblem.hpp:151
const ThermalConductionLawParams & thermalConductionLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:920
@ numComponents
Definition: FlowProblem.hpp:118
AquiferModel aquiferModel_
Definition: FlowProblem.hpp:1779
GlobalEqVector drift_
Definition: FlowProblem.hpp:1776
bool updateMinPressure_()
Definition: FlowProblem.hpp:1324
std::function< std::vector< double >(const FieldPropsManager &, const std::string &)> fieldPropDoubleOnLeafAssigner_()
Definition: FlowProblem.hpp:1357
Scalar transmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:590
@ enableExtbo
Definition: FlowProblem.hpp:128
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: FlowProblem.hpp:148
void updateReferencePorosity_()
Definition: FlowProblem.hpp:1438
Scalar thermalHalfTransmissibility(const unsigned globalSpaceIdxIn, const unsigned globalSpaceIdxOut) const
Definition: FlowProblem.hpp:610
@ dim
Definition: FlowProblem.hpp:112
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Definition: FlowProblem.hpp:1734
BCData< int > bcindex_
Definition: FlowProblem.hpp:1824
GetPropType< TypeTag, Properties::TracerModel > TracerModel
Definition: FlowProblem.hpp:168
@ enableMICP
Definition: FlowProblem.hpp:130
Scalar rockCompressibility(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:712
bool updateMaxWaterSaturation_()
Definition: FlowProblem.hpp:1293
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:935
@ numPhases
Definition: FlowProblem.hpp:117
bool updateMinPressure_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1339
std::string name() const
The problem name.
Definition: FlowProblem.hpp:884
int episodeIndex() const
Definition: FlowProblem.hpp:294
LhsEval rockCompTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1129
void endIteration()
Called by the simulator after each Newton-Raphson iteration.
Definition: FlowProblem.hpp:405
GetPropType< TypeTag, Properties::Indices > Indices
Definition: FlowProblem.hpp:109
@ numEq
Definition: FlowProblem.hpp:116
FlowProblem(Simulator &simulator)
Definition: FlowProblem.hpp:217
virtual void handleSolventBC(const BCProp::BCFace &, RateVector &) const =0
@ oilCompIdx
Definition: FlowProblem.hpp:144
GetPropType< TypeTag, Properties::GlobalEqVector > GlobalEqVector
Definition: FlowProblem.hpp:106
@ enablePolymer
Definition: FlowProblem.hpp:131
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:656
void source(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:1005
virtual void readEquilInitialCondition_()=0
Scalar nextTimeStepSize() const
Propose the size of the next time step to the simulator.
Definition: FlowProblem.hpp:1058
LhsEval computeRockCompTransMultiplier_(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the transmissibility multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1719
std::pair< BCType, RateVector > boundaryCondition(const unsigned int globalSpaceIdx, const int directionId) const
Definition: FlowProblem.hpp:1163
typename EclMaterialLawManager::MaterialLawParams MaterialLawParams
Definition: FlowProblem.hpp:154
virtual ~FlowProblem()=default
PffGridVector< GridView, Stencil, PffDofData_, DofMapper > pffDofData_
Definition: FlowProblem.hpp:1781
const ThermalConductionLawParams & thermalConductionLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:788
Scalar transmissibility(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:532
@ enableConvectiveMixing
Definition: FlowProblem.hpp:122
@ enableThermalFluxBoundaries
Definition: FlowProblem.hpp:135
TracerModel tracerModel_
Definition: FlowProblem.hpp:1782
virtual void handlePolymerBC(const BCProp::BCFace &, RateVector &) const =0
const TracerModel & tracerModel() const
Definition: FlowProblem.hpp:660
WellModel wellModel_
Definition: FlowProblem.hpp:1778
virtual void beginEpisode()
Called by the simulator before an episode begins.
Definition: FlowProblem.hpp:302
const SolidEnergyLawParams & solidEnergyLawParams(unsigned globalSpaceIdx, unsigned) const
Definition: FlowProblem.hpp:914
@ enableSaltPrecipitation
Definition: FlowProblem.hpp:133
Scalar getRockCompTransMultVal(std::size_t dofIdx) const
Definition: FlowProblem.hpp:1550
virtual void beginTimeStep()
Called by the simulator before each time integration.
Definition: FlowProblem.hpp:361
@ enableExperiments
Definition: FlowProblem.hpp:127
Scalar dofCenterDepth(unsigned globalSpaceIdx) const
Direct indexed acces to the depth of an degree of freedom [m].
Definition: FlowProblem.hpp:703
@ enableFoam
Definition: FlowProblem.hpp:129
typename GetProp< TypeTag, Properties::MaterialLaw >::EclMaterialLawManager EclMaterialLawManager
Definition: FlowProblem.hpp:152
const MaterialLawParams & materialLawParams(unsigned globalDofIdx, FaceDir::DirEnum facedir) const
Definition: FlowProblem.hpp:765
std::shared_ptr< const EclMaterialLawManager > materialLawManager() const
Returns the ECL material law manager.
Definition: FlowProblem.hpp:800
unsigned plmixnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:869
const MaterialLawParams & materialLawParams(unsigned globalDofIdx) const
Definition: FlowProblem.hpp:760
Scalar temperature(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:891
Scalar thermalHalfTransmissibilityBoundary(const Context &elemCtx, unsigned boundaryFaceIdx) const
Definition: FlowProblem.hpp:646
unsigned miscnumRegionIndex(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the index of the relevant region for thermodynmic properties.
Definition: FlowProblem.hpp:861
void updateRockCompTransMultVal_()
Definition: FlowProblem.hpp:1701
std::shared_ptr< EclThermalLawManager > thermalLawManager_
Definition: FlowProblem.hpp:1774
Scalar limitNextTimeStepSize_(Scalar dtNext) const
Definition: FlowProblem.hpp:1643
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:511
bool updateHysteresis_()
Definition: FlowProblem.hpp:1526
void readThermalParameters_()
Definition: FlowProblem.hpp:1421
void serializeOp(Serializer &serializer)
Definition: FlowProblem.hpp:1218
@ gasPhaseIdx
Definition: FlowProblem.hpp:137
Scalar maxPolymerAdsorption(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Returns the max polymer adsorption value.
Definition: FlowProblem.hpp:878
Scalar thermalHalfTransmissibilityIn(const Context &context, unsigned faceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:620
void setMaxOilSaturation(unsigned globalDofIdx, Scalar value)
Sets an element's maximum oil phase saturation observed during the simulation.
Definition: FlowProblem.hpp:952
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:775
typename EclThermalLawManager::ThermalConductionLawParams ThermalConductionLawParams
Definition: FlowProblem.hpp:156
Scalar dispersivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:570
AquiferModel & mutableAquiferModel()
Definition: FlowProblem.hpp:1046
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: FlowProblem.hpp:160
@ enableDiffusion
Definition: FlowProblem.hpp:123
Scalar temperature(unsigned globalDofIdx, unsigned) const
Definition: FlowProblem.hpp:903
@ enableBrine
Definition: FlowProblem.hpp:121
GetPropType< TypeTag, Properties::AquiferModel > AquiferModel
Definition: FlowProblem.hpp:162
TemperatureModel temperatureModel_
Definition: FlowProblem.hpp:1783
std::shared_ptr< EclMaterialLawManager > materialLawManager_
Definition: FlowProblem.hpp:1773
WellModel & wellModel()
Definition: FlowProblem.hpp:1040
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:1242
LhsEval rockCompTransMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx, Callback &obtain) const
Definition: FlowProblem.hpp:1143
bool updateMaxOilSaturation_()
Definition: FlowProblem.hpp:1261
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition: FlowProblem.hpp:182
void updatePffDofData_()
Definition: FlowProblem.hpp:1569
virtual void endEpisode()
Called by the simulator after the end of an episode.
Definition: FlowProblem.hpp:472
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: FlowProblem.hpp:147
bool nonTrivialBoundaryConditions_
Definition: FlowProblem.hpp:1825
@ waterCompIdx
Definition: FlowProblem.hpp:145
GetPropType< TypeTag, Properties::Problem > Implementation
Definition: FlowProblem.hpp:100
void readBoundaryConditions_()
Definition: FlowProblem.hpp:1600
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:804
Vanguard::TransmissibilityType transmissibilities_
Definition: FlowProblem.hpp:1771
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:996
virtual void addToSourceDense(RateVector &rate, unsigned globalDofIdx, unsigned timeIdx) const =0
virtual void endTimeStep()
Called by the simulator after each time integration.
Definition: FlowProblem.hpp:415
TemperatureModel & temperatureModel()
Definition: FlowProblem.hpp:666
Utility::CopyablePtr< DirectionalMobility< TypeTag > > DirectionalMobilityPtr
Definition: FlowProblem.hpp:169
virtual void readInitialCondition_()
Definition: FlowProblem.hpp:1497
virtual void initialSolutionApplied()
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition: FlowProblem.hpp:963
const MaterialLawParams & materialLawParams(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:753
Scalar thermalTransmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:577
static constexpr EnergyModules energyModuleType
Definition: FlowProblem.hpp:125
Scalar diffusivity(const Context &context, unsigned fromDofLocalIdx, unsigned toDofLocalIdx) const
Definition: FlowProblem.hpp:552
@ enableFullyImplicitThermal
Definition: FlowProblem.hpp:126
Scalar rockReferencePressure(unsigned globalSpaceIdx) const
Definition: FlowProblem.hpp:731
void deserialize(Restarter &res)
This method restores the complete state of the problem and its sub-objects from disk.
Definition: FlowProblem.hpp:268
std::shared_ptr< EclMaterialLawManager > materialLawManager()
Definition: FlowProblem.hpp:837
@ waterPhaseIdx
Definition: FlowProblem.hpp:139
bool updateMaxWaterSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1310
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:691
bool updateMaxOilSaturation_(unsigned compressedDofIdx, const IntensiveQuantities &iq)
Definition: FlowProblem.hpp:1279
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:1043
@ dimWorld
Definition: FlowProblem.hpp:113
@ enableDispersion
Definition: FlowProblem.hpp:124
int refPressurePhaseIdx_() const
Definition: FlowProblem.hpp:1689
LhsEval rockCompPoroMultiplier(const IntensiveQuantities &intQuants, unsigned elementIdx) const
Calculate the porosity multiplier due to water induced rock compaction.
Definition: FlowProblem.hpp:1085
Scalar transmissibilityBoundary(const unsigned globalSpaceIdx, const unsigned boundaryFaceIdx) const
Direct access to a boundary transmissibility.
Definition: FlowProblem.hpp:600
TracerModel & tracerModel()
Definition: FlowProblem.hpp:663
void readMaterialParameters_()
Definition: FlowProblem.hpp:1381
void serialize(Restarter &res)
This method writes the complete state of the problem and its subobjects to disk.
Definition: FlowProblem.hpp:287
MathToolbox< Evaluation > Toolbox
Definition: FlowProblem.hpp:164
Scalar rockReferencePressure(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition: FlowProblem.hpp:722
Scalar diffusivity(const unsigned globalCellIn, const unsigned globalCellOut) const
Definition: FlowProblem.hpp:563
void updateRockFraction_()
Definition: FlowProblem.hpp:1463
void prefetch(const Element &elem) const
Definition: FlowProblem.hpp:253
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:1787
const std::vector< T > & operator()(FaceDir::DirEnum dir) const
Definition: FlowProblem.hpp:1796
void resize(std::size_t size, T defVal)
Definition: FlowProblem.hpp:1790
std::vector< T > & operator()(FaceDir::DirEnum dir)
Definition: FlowProblem.hpp:1808
std::array< std::vector< T >, 6 > data
Definition: FlowProblem.hpp:1788
Definition: FlowProblem.hpp:1560
ConditionalStorage< enableFullyImplicitThermal, Scalar > thermalHalfTransOut
Definition: FlowProblem.hpp:1562
ConditionalStorage< enableFullyImplicitThermal, Scalar > thermalHalfTransIn
Definition: FlowProblem.hpp:1561
ConditionalStorage< enableDiffusion, Scalar > diffusivity
Definition: FlowProblem.hpp:1563
ConditionalStorage< enableDispersion, Scalar > dispersivity
Definition: FlowProblem.hpp:1564
Scalar transmissibility
Definition: FlowProblem.hpp:1565