opm-simulators
Loading...
Searching...
No Matches
BlackoilWellModel.hpp
1/*
2 Copyright 2016 SINTEF ICT, Applied Mathematics.
3 Copyright 2016 - 2017 Statoil ASA.
4 Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
5 Copyright 2016 - 2018 IRIS AS
6
7 This file is part of the Open Porous Media project (OPM).
8
9 OPM is free software: you can redistribute it and/or modify
10 it under the terms of the GNU General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 OPM is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18
19 You should have received a copy of the GNU General Public License
20 along with OPM. If not, see <http://www.gnu.org/licenses/>.
21*/
22
23#ifndef OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
24#define OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
25
26#if HAVE_MPI
27#define RESERVOIR_COUPLING_ENABLED
28#endif
29#ifdef RESERVOIR_COUPLING_ENABLED
30#include <opm/simulators/flow/ReservoirCouplingMaster.hpp>
31#include <opm/simulators/flow/ReservoirCouplingSlave.hpp>
32#endif
33
34#include <dune/common/fmatrix.hh>
35#include <dune/istl/bcrsmatrix.hh>
36#include <dune/istl/matrixmatrix.hh>
37
38#include <opm/common/OpmLog/OpmLog.hpp>
39
40#include <opm/input/eclipse/Schedule/Group/Group.hpp>
41#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
42#include <opm/input/eclipse/Schedule/Schedule.hpp>
43#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
44
45#include <opm/material/densead/Math.hpp>
46
47#include <opm/simulators/flow/countGlobalCells.hpp>
49
50#include <opm/simulators/linalg/matrixblock.hh>
51
52#include <opm/simulators/timestepping/SimulatorReport.hpp>
53#include <opm/simulators/timestepping/gatherConvergenceReport.hpp>
54
55#include <opm/simulators/utils/DeferredLogger.hpp>
56
57#include <opm/simulators/wells/BlackoilWellModelGasLift.hpp>
58#include <opm/simulators/wells/BlackoilWellModelGeneric.hpp>
59#include <opm/simulators/wells/BlackoilWellModelGuideRates.hpp>
60#include <opm/simulators/wells/GasLiftGroupInfo.hpp>
61#include <opm/simulators/wells/GasLiftSingleWell.hpp>
62#include <opm/simulators/wells/GasLiftSingleWellGeneric.hpp>
63#include <opm/simulators/wells/GasLiftWellState.hpp>
64#include <opm/simulators/wells/GuideRateHandler.hpp>
65#include <opm/simulators/wells/MultisegmentWell.hpp>
66#include <opm/simulators/wells/ParallelWBPCalculation.hpp>
67#include <opm/simulators/wells/ParallelWellInfo.hpp>
68#include <opm/simulators/wells/PerforationData.hpp>
71#include <opm/simulators/wells/StandardWell.hpp>
72#include <opm/simulators/wells/VFPInjProperties.hpp>
73#include <opm/simulators/wells/VFPProdProperties.hpp>
74#include <opm/simulators/wells/WellConnectionAuxiliaryModule.hpp>
75#include <opm/simulators/wells/WellGroupHelpers.hpp>
76#include <opm/simulators/wells/WellInterface.hpp>
77#include <opm/simulators/wells/WellProdIndexCalculator.hpp>
78#include <opm/simulators/wells/WellState.hpp>
79#include <opm/simulators/wells/WGState.hpp>
80
81#include <cstddef>
82#include <map>
83#include <memory>
84#include <string>
85#include <tuple>
86#include <vector>
87
88namespace Opm {
89
90template<class Scalar> class BlackoilWellModelNldd;
91template<class T> class SparseTable;
92
93#if COMPILE_GPU_BRIDGE
94template<class Scalar> class WellContributions;
95#endif
96
98 template<typename TypeTag>
99 class BlackoilWellModel : public WellConnectionAuxiliaryModule<TypeTag, BlackoilWellModel<TypeTag>>
100 , public BlackoilWellModelGeneric<GetPropType<TypeTag, Properties::Scalar>,
101 typename GetPropType<TypeTag, Properties::FluidSystem>::IndexTraitsType>
102 {
103 public:
104 // --------- Types ---------
115 using ModelParameters = BlackoilModelParameters<Scalar>;
116
117 using WellConnectionModule = WellConnectionAuxiliaryModule<TypeTag, BlackoilWellModel<TypeTag>>;
118
119 using IndexTraits = typename FluidSystem::IndexTraitsType;
120
121 constexpr static std::size_t pressureVarIndex = GetPropType<TypeTag, Properties::Indices>::pressureSwitchIdx;
122
123 static const int numEq = Indices::numEq;
124 static const int solventSaturationIdx = Indices::solventSaturationIdx;
125 static constexpr bool has_solvent_ = getPropValue<TypeTag, Properties::EnableSolvent>();
126 static constexpr bool has_polymer_ = getPropValue<TypeTag, Properties::EnablePolymer>();
127 static constexpr bool has_energy_ = getPropValue<TypeTag, Properties::EnableEnergy>();
128 static constexpr bool has_micp_ = Indices::enableMICP;
129
130 // TODO: where we should put these types, WellInterface or Well Model?
131 // or there is some other strategy, like TypeTag
132 using VectorBlockType = Dune::FieldVector<Scalar, numEq>;
133 using BVector = Dune::BlockVector<VectorBlockType>;
134
135 using PolymerModule = BlackOilPolymerModule<TypeTag>;
136 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag>;
137
138 // For the conversion between the surface volume rate and reservoir voidage rate
139 using RateConverterType = RateConverter::
140 SurfaceToReservoirVoidage<FluidSystem, std::vector<int> >;
141
142 // For computing average pressured used by gpmaint
143 using AverageRegionalPressureType = RegionAverageCalculator::
144 AverageRegionalPressure<FluidSystem, std::vector<int> >;
145
146 using WellGroupHelpersType = WellGroupHelpers<Scalar, IndexTraits>;
147
148 explicit BlackoilWellModel(Simulator& simulator);
149
150 void init();
151 void initWellContainer(const int reportStepIdx) override;
152
153 void beginEpisode()
154 {
155 OPM_TIMEBLOCK(beginEpsiode);
156 beginReportStep(simulator_.episodeIndex());
157 }
158
159 void beginTimeStep();
160
161 void beginIteration()
162 {
163 OPM_TIMEBLOCK(beginIteration);
164 assemble(simulator_.model().newtonMethod().numIterations(),
165 simulator_.timeStepSize());
166 }
167
168 void endIteration()
169 { }
170
171 void endTimeStep()
172 {
173 OPM_TIMEBLOCK(endTimeStep);
174 timeStepSucceeded(simulator_.time(), simulator_.timeStepSize());
175 }
176
177 void endEpisode()
178 {
179 endReportStep();
180 }
181
182 void computeTotalRatesForDof(RateVector& rate,
183 unsigned globalIdx) const;
184
185 template <class Context>
186 void computeTotalRatesForDof(RateVector& rate,
187 const Context& context,
188 unsigned spaceIdx,
189 unsigned timeIdx) const;
190
191
192 using WellInterfacePtr = std::shared_ptr<WellInterface<TypeTag> >;
193
194 using BlackoilWellModelGeneric<Scalar, IndexTraits>::initFromRestartFile;
195 void initFromRestartFile(const RestartValue& restartValues)
196 {
197 initFromRestartFile(restartValues,
198 this->simulator_.vanguard().transferWTestState(),
199 grid().size(0),
200 param_.use_multisegment_well_,
201 this->simulator_.vanguard().enableDistributedWells());
202 }
203
204 using BlackoilWellModelGeneric<Scalar, IndexTraits>::prepareDeserialize;
205 void prepareDeserialize(const int report_step)
206 {
207 prepareDeserialize(report_step, grid().size(0),
208 param_.use_multisegment_well_,
209 this->simulator_.vanguard().enableDistributedWells());
210 }
211
212 data::Wells wellData() const
213 {
214 auto wsrpt = this->wellState()
215 .report(this->simulator_.vanguard().globalCell().data(),
216 [this](const int well_index)
217 {
218 return this->wasDynamicallyShutThisTimeStep(well_index);
219 });
220
222 .assignWellGuideRates(wsrpt, this->reportStepIndex());
223
224 this->assignWellTracerRates(wsrpt);
225
226 if (const auto& rspec = eclState().runspec();
227 rspec.co2Storage() || rspec.h2Storage())
228 {
229 // The gas reference density (surface condition) is the
230 // same for all PVT regions in CO2STORE/H2STORE runs so,
231 // for simplicity, we use region zero (0) here.
232
233 this->assignMassGasRate(wsrpt, FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, 0));
234 }
235
236 this->assignWellTargets(wsrpt);
237
238 this->assignDynamicWellStatus(wsrpt);
239
240 // Assigning (a subset of the) property values in shut
241 // connections should be the last step of wellData().
242 this->assignShutConnections(wsrpt, this->reportStepIndex());
243
244 return wsrpt;
245 }
246
247 data::WellBlockAveragePressures wellBlockAveragePressures() const
248 {
249 return this->wbp_.computeWellBlockAveragePressures(this->gravity_);
250 }
251
252#if COMPILE_GPU_BRIDGE
253 // accumulate the contributions of all Wells in the WellContributions object
254 void getWellContributions(WellContributions<Scalar>& x) const;
255#endif
256
257 // Check if well equations is converged.
258 ConvergenceReport getWellConvergence(const std::vector<Scalar>& B_avg, const bool checkWellGroupControlsAndNetwork = false) const;
259
260 const SimulatorReportSingle& lastReport() const;
261
262 void addWellContributions(SparseMatrixAdapter& jacobian) const;
263
264 // add source from wells to the reservoir matrix
265 void addReservoirSourceTerms(GlobalEqVector& residual,
266 const std::vector<typename SparseMatrixAdapter::MatrixBlock*>& diagMatAddress) const;
267
268 // called at the beginning of a report step
269 void beginReportStep(const int time_step);
270
271 // it should be able to go to prepareTimeStep(), however, the updateWellControls()
272 // makes it a little more difficult. unless we introduce if (iterationIdx != 0) to avoid doing the above function
273 // twice at the beginning of the time step
276 void calculateExplicitQuantities(DeferredLogger& deferred_logger) const;
277 // some preparation work, mostly related to group control and RESV,
278 // at the beginning of each time step (Not report step)
279 void prepareTimeStep(DeferredLogger& deferred_logger);
280
281 bool
282 updateWellControls(DeferredLogger& deferred_logger);
283
284 std::tuple<bool, Scalar>
285 updateNetworks(const bool mandatory_network_balance, DeferredLogger& deferred_logger, const bool relax_network_tolerance = false);
286
287
288 void updateAndCommunicate(const int reportStepIdx,
289 const int iterationIdx,
290 DeferredLogger& deferred_logger);
291
292 bool updateGroupControls(const Group& group,
293 DeferredLogger& deferred_logger,
294 const int reportStepIdx,
295 const int iterationIdx);
296
297 WellInterfacePtr getWell(const std::string& well_name) const;
298
299 using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
300
301 void addWellPressureEquations(PressureMatrix& jacobian,
302 const BVector& weights,
303 const bool use_well_weights) const;
304 void addWellPressureEquationsStruct(PressureMatrix& jacobian) const;
305 void addWellPressureEquationsDomain(PressureMatrix& jacobian,
306 const BVector& weights,
307 const bool use_well_weights,
308 const int domainIndex) const
309 {
310 if (!nldd_) {
311 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
312 }
313 return nldd_->addWellPressureEquations(jacobian,
314 weights,
315 use_well_weights,
316 domainIndex);
317 }
318
320 const std::vector<WellInterfacePtr>& localNonshutWells() const
321 {
322 return well_container_;
323 }
324
325 const SparseTable<int>& well_local_cells() const
326 {
327 if (!nldd_) {
328 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
329 }
330 return nldd_->well_local_cells();
331 }
332
333 const std::map<std::string, int>& well_domain() const
334 {
335 if (!nldd_) {
336 OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
337 }
338
339 return nldd_->well_domain();
340 }
341
342 auto begin() const { return well_container_.begin(); }
343 auto end() const { return well_container_.end(); }
344 bool empty() const { return well_container_.empty(); }
345
346 bool addMatrixContributions() const
347 { return param_.matrix_add_well_contributions_; }
348
349 int numStrictIterations() const
350 { return param_.strict_outer_iter_wells_; }
351
352 int compressedIndexForInterior(int cartesian_cell_idx) const override
353 {
354 return simulator_.vanguard().compressedIndexForInterior(cartesian_cell_idx);
355 }
356
357 int compressedIndexForInteriorLGR(const std::string& lgr_tag, const Connection& conn) const override
358 {
359 return simulator_.vanguard().compressedIndexForInteriorLGR(lgr_tag, conn);
360 }
361
362 // using the solution x to recover the solution xw for wells and applying
363 // xw to update Well State
364 void recoverWellSolutionAndUpdateWellState(const BVector& x);
365
366 // using the solution x to recover the solution xw for wells and applying
367 // xw to update Well State
368 void recoverWellSolutionAndUpdateWellStateDomain(const BVector& x,
369 const int domainIdx);
370
371 const Grid& grid() const
372 { return simulator_.vanguard().grid(); }
373
374 const Simulator& simulator() const
375 { return simulator_; }
376
377 void setNlddAdapter(BlackoilWellModelNldd<TypeTag>* mod)
378 { nldd_ = mod; }
379
380#ifdef RESERVOIR_COUPLING_ENABLED
381 ReservoirCouplingMaster& reservoirCouplingMaster() {
382 return *(this->simulator_.reservoirCouplingMaster());
383 }
384 ReservoirCouplingSlave& reservoirCouplingSlave() {
385 return *(this->simulator_.reservoirCouplingSlave());
386 }
387 bool isReservoirCouplingMaster() const {
388 return this->simulator_.reservoirCouplingMaster() != nullptr;
389 }
390 bool isReservoirCouplingSlave() const {
391 return this->simulator_.reservoirCouplingSlave() != nullptr;
392 }
393 void setReservoirCouplingMaster(ReservoirCouplingMaster *master)
394 {
395 this->guide_rate_handler_.setReservoirCouplingMaster(master);
396 }
397 void setReservoirCouplingSlave(ReservoirCouplingSlave *slave)
398 {
399 this->guide_rate_handler_.setReservoirCouplingSlave(slave);
400 }
401 #endif
402 protected:
403 Simulator& simulator_;
404
405 // a vector of all the wells.
406 std::vector<WellInterfacePtr> well_container_{};
407
408 std::vector<bool> is_cell_perforated_{};
409
410 void initializeWellState(const int timeStepIdx);
411
412 // create the well container
413 void createWellContainer(const int report_step) override;
414
415 WellInterfacePtr
416 createWellPointer(const int wellID,
417 const int report_step) const;
418
419 template <typename WellType>
420 std::unique_ptr<WellType>
421 createTypedWellPointer(const int wellID,
422 const int time_step) const;
423
424 WellInterfacePtr createWellForWellTest(const std::string& well_name,
425 const int report_step,
426 DeferredLogger& deferred_logger) const;
427
428 const ModelParameters param_;
429 std::size_t global_num_cells_{};
430 // the number of the cells in the local grid
431 std::size_t local_num_cells_{};
432 Scalar gravity_{};
433 std::vector<Scalar> depth_{};
434 bool alternative_well_rate_init_{};
435 std::map<std::string, Scalar> well_group_thp_calc_;
436 std::unique_ptr<RateConverterType> rateConverter_{};
437 std::map<std::string, std::unique_ptr<AverageRegionalPressureType>> regionalAveragePressureCalculator_{};
438
439 SimulatorReportSingle last_report_{};
440 GuideRateHandler<Scalar, IndexTraits> guide_rate_handler_{};
441
442 // A flag to tell the convergence report whether we need to take another newton step
443 bool network_needs_more_balancing_force_another_newton_iteration_{false};
444
445 // Pre-step network solve at static reservoir conditions (group and well states might be updated)
446 void doPreStepNetworkRebalance(DeferredLogger& deferred_logger);
447
448 std::vector<Scalar> B_avg_{};
449
450 const EquilGrid& equilGrid() const
451 { return simulator_.vanguard().equilGrid(); }
452
453 const EclipseState& eclState() const
454 { return simulator_.vanguard().eclState(); }
455
456 // compute the well fluxes and assemble them in to the reservoir equations as source terms
457 // and in the well equations.
458 void assemble(const int iterationIdx,
459 const double dt);
460
461 // well controls and network pressures affect each other and are solved in an iterative manner.
462 // the function handles one iteration of updating well controls and network pressures.
463 // it is possible to decouple the update of well controls and network pressures further.
464 // the returned two booleans are {continue_due_to_network, well_group_control_changed}, respectively
465 std::tuple<bool, bool, Scalar> updateWellControlsAndNetworkIteration(const bool mandatory_network_balance,
466 const bool relax_network_tolerance,
467 const bool optimize_gas_lift,
468 const double dt,
469 DeferredLogger& local_deferredLogger);
470
471 bool updateWellControlsAndNetwork(const bool mandatory_network_balance,
472 const double dt,
473 DeferredLogger& local_deferredLogger);
474
475 bool computeWellGroupThp(const double dt, DeferredLogger& local_deferredLogger);
476
485 void initializeLocalWellStructure(const int reportStepIdx,
486 const bool enableWellPIScaling);
487
491 void initializeGroupStructure(const int reportStepIdx);
492
493 // called at the end of a time step
494 void timeStepSucceeded(const double simulationTime, const double dt);
495
496 // called at the end of a report step
497 void endReportStep();
498
499 // setting the well_solutions_ based on well_state.
500 void updatePrimaryVariables(DeferredLogger& deferred_logger);
501
502 void updateAverageFormationFactor();
503
504 void computePotentials(const std::size_t widx,
505 const WellState<Scalar, IndexTraits>& well_state_copy,
506 std::string& exc_msg,
507 ExceptionType::ExcEnum& exc_type,
508 DeferredLogger& deferred_logger) override;
509
510 const std::vector<Scalar>& wellPerfEfficiencyFactors() const;
511
512 void calculateProductivityIndexValuesShutWells(const int reportStepIdx, DeferredLogger& deferred_logger) override;
513 void calculateProductivityIndexValues(DeferredLogger& deferred_logger) override;
514 void calculateProductivityIndexValues(const WellInterface<TypeTag>* wellPtr,
515 DeferredLogger& deferred_logger);
516
517 // The number of conservation quantities.
518 int numConservationQuantities() const;
519
520 int reportStepIndex() const;
521
522 void assembleWellEq(const double dt, DeferredLogger& deferred_logger);
523
524 void prepareWellsBeforeAssembling(const double dt, DeferredLogger& deferred_logger);
525
526 // TODO: finding a better naming
527 void assembleWellEqWithoutIteration(const double dt, DeferredLogger& deferred_logger);
528
529 void extractLegacyCellPvtRegionIndex_();
530
531 void extractLegacyDepth_();
532
534 void updateWellTestState(const double simulationTime, WellTestState& wellTestState);
535
536 void wellTesting(const int timeStepIdx, const double simulationTime, DeferredLogger& deferred_logger);
537
538 void calcResvCoeff(const int fipnum,
539 const int pvtreg,
540 const std::vector<Scalar>& production_rates,
541 std::vector<Scalar>& resv_coeff) override;
542
543 void calcInjResvCoeff(const int fipnum,
544 const int pvtreg,
545 std::vector<Scalar>& resv_coeff) override;
546
547 void computeWellTemperature();
548
549 private:
550 BlackoilWellModelGasLift<TypeTag> gaslift_;
551 BlackoilWellModelNldd<TypeTag>* nldd_ = nullptr;
552
553 // These members are used to avoid reallocation in specific functions
554 // instead of using local variables.
555 // Their state is not relevant between function calls, so they can
556 // (and must) be mutable, as the functions using them are const.
557 mutable BVector x_local_;
558
559 void assignWellTracerRates(data::Wells& wsrpt) const;
560 };
561
562} // namespace Opm
563
564#include "BlackoilWellModel_impl.hpp"
565
566#endif // OPM_BLACKOILWELLMODEL_HEADER_INCLUDED
Helper class for grid instantiation of ECL file-format using problems.
Facility for converting component rates at surface conditions to phase (voidage) rates at reservoir c...
Facility for converting component rates at surface conditions to phase (voidage) rates at reservoir c...
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition blackoilbioeffectsmodules.hh:93
Contains the high level supplements required to extend the black oil model by polymer.
Definition blackoilpolymermodules.hh:64
void assignShutConnections(data::Wells &wsrpt, const int reportStepIndex) const
Definition BlackoilWellModelGeneric.cpp:1154
Class for handling the guide rates in the blackoil well model.
Definition BlackoilWellModelGuideRates.hpp:46
void assignWellGuideRates(data::Wells &wsrpt, const int reportStepIdx) const
Assign well guide rates.
Definition BlackoilWellModelGuideRates.cpp:481
Class for handling the blackoil well model in a NLDD solver.
Definition BlackoilWellModelNldd.hpp:80
void initializeGroupStructure(const int reportStepIdx)
Initialize group control modes/constraints and group solution state.
Definition BlackoilWellModel_impl.hpp:291
const std::vector< WellInterfacePtr > & localNonshutWells() const
Get list of local nonshut wells.
Definition BlackoilWellModel.hpp:320
void calculateExplicitQuantities(DeferredLogger &deferred_logger) const
Calculating the explict quantities used in the well calculation.
Definition BlackoilWellModel_impl.hpp:1767
void initializeLocalWellStructure(const int reportStepIdx, const bool enableWellPIScaling)
Update rank's notion of intersecting wells and their associate solution variables.
Definition BlackoilWellModel_impl.hpp:246
void updateWellTestState(const double simulationTime, WellTestState &wellTestState)
upate the wellTestState related to economic limits
Definition BlackoilWellModel_impl.hpp:2001
int compressedIndexForInterior(int cartesian_cell_idx) const override
get compressed index for interior cells (-1, otherwise
Definition BlackoilWellModel.hpp:352
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition ConvergenceReport.hpp:38
Definition DeferredLogger.hpp:57
Definition BlackoilWellModel.hpp:91
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51
Definition WellGroupHelpers.hpp:51
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
Solver parameters for the BlackoilModel.
Definition BlackoilModelParameters.hpp:180
A struct for returning timing data from a simulator to its caller.
Definition SimulatorReport.hpp:34