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