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