WellInterface_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4 Copyright 2018 IRIS
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 3 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
22// Improve IDE experience
23#ifndef OPM_WELLINTERFACE_HEADER_INCLUDED
24#include <config.h>
25#define OPM_WELLINTERFACE_IMPL_HEADER_INCLUDED
27#endif
28
29#include <opm/common/Exceptions.hpp>
30
31#include <opm/input/eclipse/Schedule/ScheduleTypes.hpp>
32#include <opm/input/eclipse/Schedule/Well/WDFAC.hpp>
33
35
40
41#include <dune/common/version.hh>
42
43#include <cassert>
44#include <cstddef>
45#include <utility>
46
47#include <fmt/format.h>
48
49namespace Opm
50{
51
52
53 template<typename TypeTag>
55 WellInterface(const Well& well,
56 const ParallelWellInfo<Scalar>& pw_info,
57 const int time_step,
58 const ModelParameters& param,
59 const RateConverterType& rate_converter,
60 const int pvtRegionIdx,
61 const int num_components,
62 const int num_phases,
63 const int index_of_well,
64 const std::vector<PerforationData<Scalar>>& perf_data)
66 pw_info,
67 time_step,
68 rate_converter,
69 pvtRegionIdx,
70 num_components,
71 num_phases,
72 index_of_well,
73 perf_data)
74 , param_(param)
75 {
77
78 if constexpr (has_solvent || has_zFraction) {
79 if (well.isInjector()) {
80 auto injectorType = this->well_ecl_.injectorType();
81 if (injectorType == InjectorType::GAS) {
82 this->wsolvent_ = this->well_ecl_.getSolventFraction();
83 }
84 }
85 }
86 }
87
88
89 template<typename TypeTag>
90 void
92 init(const PhaseUsage* phase_usage_arg,
93 const std::vector<Scalar>& /* depth_arg */,
94 const Scalar gravity_arg,
95 const int /* num_cells */,
96 const std::vector<Scalar>& B_avg,
97 const bool changed_to_open_this_step)
98 {
99 this->phase_usage_ = phase_usage_arg;
100 this->gravity_ = gravity_arg;
101 B_avg_ = B_avg;
102 this->changed_to_open_this_step_ = changed_to_open_this_step;
103 }
104
105
106
107
108 template<typename TypeTag>
111 wpolymer() const
112 {
113 if constexpr (has_polymer) {
114 return this->wpolymer_();
115 }
116
117 return 0.0;
118 }
119
120
121
122
123
124 template<typename TypeTag>
127 wfoam() const
128 {
129 if constexpr (has_foam) {
130 return this->wfoam_();
131 }
132
133 return 0.0;
134 }
135
136
137
138 template<typename TypeTag>
141 wsalt() const
142 {
143 if constexpr (has_brine) {
144 return this->wsalt_();
145 }
146
147 return 0.0;
148 }
149
150 template<typename TypeTag>
153 wmicrobes() const
154 {
155 if constexpr (has_micp) {
156 return this->wmicrobes_();
157 }
158
159 return 0.0;
160 }
161
162 template<typename TypeTag>
165 woxygen() const
166 {
167 if constexpr (has_micp) {
168 return this->woxygen_();
169 }
170
171 return 0.0;
172 }
173
174 // The urea injection concentration is scaled down by a factor of 10, since its value
175 // can be much bigger than 1 (not doing this slows the simulations). The
176 // corresponding values are scaled accordingly in blackoilmicpmodules.hh when computing
177 // the reactions and also when writing the output files (vtk and eclipse format, i.e.,
178 // vtkblackoilmicpmodule.hh and ecloutputblackoilmodel.hh respectively).
179
180 template<typename TypeTag>
183 wurea() const
184 {
185 if constexpr (has_micp) {
186 return this->wurea_();
187 }
188
189 return 0.0;
190 }
191
192 template<typename TypeTag>
193 bool
195 updateWellControl(const Simulator& simulator,
196 const IndividualOrGroup iog,
197 WellState<Scalar>& well_state,
198 const GroupState<Scalar>& group_state,
199 DeferredLogger& deferred_logger) /* const */
200 {
201 if (stoppedOrZeroRateTarget(simulator, well_state, deferred_logger)) {
202 return false;
203 }
204
205 const auto& summaryState = simulator.vanguard().summaryState();
206 const auto& schedule = simulator.vanguard().schedule();
207 const auto& well = this->well_ecl_;
208 auto& ws = well_state.well(this->index_of_well_);
209 std::string from;
210 if (well.isInjector()) {
211 from = WellInjectorCMode2String(ws.injection_cmode);
212 } else {
213 from = WellProducerCMode2String(ws.production_cmode);
214 }
215 bool oscillating = std::count(this->well_control_log_.begin(), this->well_control_log_.end(), from) >= param_.max_number_of_well_switches_;
216
217 if (oscillating) {
218 // only output frist time
219 bool output = std::count(this->well_control_log_.begin(), this->well_control_log_.end(), from) == param_.max_number_of_well_switches_;
220 if (output) {
221 std::ostringstream ss;
222 ss << " The control mode for well " << this->name()
223 << " is oscillating\n"
224 << " We don't allow for more than "
225 << param_.max_number_of_well_switches_
226 << " switches. The control is kept at " << from;
227 deferred_logger.info(ss.str());
228 // add one more to avoid outputting the same info again
229 this->well_control_log_.push_back(from);
230 }
231 return false;
232 }
233 bool changed = false;
234 if (iog == IndividualOrGroup::Individual) {
235 changed = this->checkIndividualConstraints(ws, summaryState, deferred_logger);
236 } else if (iog == IndividualOrGroup::Group) {
237 changed = this->checkGroupConstraints(well_state, group_state, schedule, summaryState, deferred_logger);
238 } else {
239 assert(iog == IndividualOrGroup::Both);
240 changed = this->checkConstraints(well_state, group_state, schedule, summaryState, deferred_logger);
241 }
242 Parallel::Communication cc = simulator.vanguard().grid().comm();
243 // checking whether control changed
244 if (changed) {
245 std::string to;
246 if (well.isInjector()) {
247 to = WellInjectorCMode2String(ws.injection_cmode);
248 } else {
249 to = WellProducerCMode2String(ws.production_cmode);
250 }
251 std::ostringstream ss;
252 ss << " Switching control mode for well " << this->name()
253 << " from " << from
254 << " to " << to;
255 if (cc.size() > 1) {
256 ss << " on rank " << cc.rank();
257 }
258 deferred_logger.debug(ss.str());
259
260 this->well_control_log_.push_back(from);
261 updateWellStateWithTarget(simulator, group_state, well_state, deferred_logger);
262 updatePrimaryVariables(simulator, well_state, deferred_logger);
263 }
264
265 return changed;
266 }
267
268 template<typename TypeTag>
269 bool
272 WellState<Scalar>& well_state,
273 const GroupState<Scalar>& group_state,
274 const Well::InjectionControls& inj_controls,
275 const Well::ProductionControls& prod_controls,
276 const Scalar wqTotal,
277 DeferredLogger& deferred_logger,
278 const bool fixed_control,
279 const bool fixed_status)
280 {
281 const auto& summary_state = simulator.vanguard().summaryState();
282 auto& ws = well_state.well(this->index_of_well_);
283 std::string from;
284 if (this->isInjector()) {
285 from = WellInjectorCMode2String(ws.injection_cmode);
286 } else {
287 from = WellProducerCMode2String(ws.production_cmode);
288 }
289 const bool oscillating = std::count(this->well_control_log_.begin(), this->well_control_log_.end(), from) >= param_.max_number_of_well_switches_;
290
291 if (oscillating || this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) || !(this->well_ecl_.getStatus() == WellStatus::OPEN)) {
292 return false;
293 }
294
295 const Scalar sgn = this->isInjector() ? 1.0 : -1.0;
296 if (!this->wellIsStopped()){
297 if (wqTotal*sgn <= 0.0 && !fixed_status){
298 this->stopWell();
299 return true;
300 } else {
301 bool changed = false;
302 if (!fixed_control) {
303 // We don't allow changing to group controls here since this may lead to inconsistencies
304 // in the group handling which in turn may result in excessive back and forth switching.
305 // If we are to allow this change, one needs to make sure all necessary information propagates
306 // properly, but for now, we simply disallow it. The commented code below is kept for future reference
307
308 // const bool hasGroupControl = this->isInjector() ? inj_controls.hasControl(Well::InjectorCMode::GRUP) :
309 // prod_controls.hasControl(Well::ProducerCMode::GRUP);
310
311 changed = this->checkIndividualConstraints(ws, summary_state, deferred_logger, inj_controls, prod_controls);
312 // if (hasGroupControl) {
313 // changed = changed || this->checkGroupConstraints(well_state, group_state, schedule, summary_state,deferred_logger);
314 // }
315
316 if (changed) {
317 const bool thp_controlled = this->isInjector() ? ws.injection_cmode == Well::InjectorCMode::THP :
318 ws.production_cmode == Well::ProducerCMode::THP;
319 if (!thp_controlled){
320 // don't call for thp since this might trigger additional local solve
321 updateWellStateWithTarget(simulator, group_state, well_state, deferred_logger);
322 } else {
323 ws.thp = this->getTHPConstraint(summary_state);
324 }
325 updatePrimaryVariables(simulator, well_state, deferred_logger);
326 }
327 }
328 return changed;
329 }
330 } else if (!fixed_status){
331 // well is stopped, check if current bhp allows reopening
332 const Scalar bhp = well_state.well(this->index_of_well_).bhp;
333 Scalar prod_limit = prod_controls.bhp_limit;
334 Scalar inj_limit = inj_controls.bhp_limit;
335 const bool has_thp = this->wellHasTHPConstraints(summary_state);
336 if (has_thp){
337 std::vector<Scalar> rates(this->num_components_);
338 if (this->isInjector()){
339 const Scalar bhp_thp = WellBhpThpCalculator(*this).
340 calculateBhpFromThp(well_state, rates,
341 this->well_ecl_,
342 summary_state,
343 this->getRefDensity(),
344 deferred_logger);
345 inj_limit = std::min(bhp_thp, static_cast<Scalar>(inj_controls.bhp_limit));
346 } else {
347 // if the well can operate, it must at least be able to produce
348 // at the lowest bhp of the bhp-curve (explicit fractions)
349 const Scalar bhp_min = WellBhpThpCalculator(*this).
350 calculateMinimumBhpFromThp(well_state,
351 this->well_ecl_,
352 summary_state,
353 this->getRefDensity());
354 prod_limit = std::max(bhp_min, static_cast<Scalar>(prod_controls.bhp_limit));
355 }
356 }
357 const Scalar bhp_diff = (this->isInjector())? inj_limit - bhp: bhp - prod_limit;
358 if (bhp_diff > 0){
359 this->openWell();
360 well_state.well(this->index_of_well_).bhp = (this->isInjector())? inj_limit : prod_limit;
361 if (has_thp) {
362 well_state.well(this->index_of_well_).thp = this->getTHPConstraint(summary_state);
363 }
364 return true;
365 } else {
366 return false;
367 }
368 } else {
369 return false;
370 }
371 }
372
373 template<typename TypeTag>
374 void
376 wellTesting(const Simulator& simulator,
377 const double simulation_time,
378 /* const */ WellState<Scalar>& well_state,
379 const GroupState<Scalar>& group_state,
380 WellTestState& well_test_state,
381 DeferredLogger& deferred_logger)
382 {
383 deferred_logger.info(" well " + this->name() + " is being tested");
384
385 WellState<Scalar> well_state_copy = well_state;
386 auto& ws = well_state_copy.well(this->indexOfWell());
387
388 updateWellStateWithTarget(simulator, group_state, well_state_copy, deferred_logger);
389 calculateExplicitQuantities(simulator, well_state_copy, deferred_logger);
390 updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
391 initPrimaryVariablesEvaluation();
392
393 if (this->isProducer()) {
394 const auto& schedule = simulator.vanguard().schedule();
395 const auto report_step = simulator.episodeIndex();
396 const auto& glo = schedule.glo(report_step);
397 if (glo.active()) {
398 gliftBeginTimeStepWellTestUpdateALQ(simulator, well_state_copy, deferred_logger);
399 }
400 }
401
402 WellTestState welltest_state_temp;
403
404 bool testWell = true;
405 // if a well is closed because all completions are closed, we need to check each completion
406 // individually. We first open all completions, then we close one by one by calling updateWellTestState
407 // untill the number of closed completions do not increase anymore.
408 while (testWell) {
409 const std::size_t original_number_closed_completions = welltest_state_temp.num_closed_completions();
410 bool converged = solveWellForTesting(simulator, well_state_copy, group_state, deferred_logger);
411 if (!converged) {
412 const auto msg = fmt::format("WTEST: Well {} is not solvable (physical)", this->name());
413 deferred_logger.debug(msg);
414 return;
415 }
416
417
418 updateWellOperability(simulator, well_state_copy, deferred_logger);
419 if ( !this->isOperableAndSolvable() ) {
420 const auto msg = fmt::format("WTEST: Well {} is not operable (physical)", this->name());
421 deferred_logger.debug(msg);
422 return;
423 }
424 std::vector<Scalar> potentials;
425 try {
426 computeWellPotentials(simulator, well_state_copy, potentials, deferred_logger);
427 } catch (const std::exception& e) {
428 const std::string msg = fmt::format("well {}: computeWellPotentials() "
429 "failed during testing for re-opening: ",
430 this->name(), e.what());
431 deferred_logger.info(msg);
432 return;
433 }
434 const int np = well_state_copy.numPhases();
435 for (int p = 0; p < np; ++p) {
436 ws.well_potentials[p] = std::max(Scalar{0.0}, potentials[p]);
437 }
438 const bool under_zero_target = this->wellUnderZeroGroupRateTarget(simulator, well_state_copy, deferred_logger);
439 this->updateWellTestState(well_state_copy.well(this->indexOfWell()),
440 simulation_time,
441 /*writeMessageToOPMLog=*/ false,
442 under_zero_target,
443 welltest_state_temp,
444 deferred_logger);
445 this->closeCompletions(welltest_state_temp);
446
447 // Stop testing if the well is closed or shut due to all completions shut
448 // Also check if number of completions has increased. If the number of closed completions do not increased
449 // we stop the testing.
450 // TODO: it can be tricky here, if the well is shut/closed due to other reasons
451 if ( welltest_state_temp.num_closed_wells() > 0 ||
452 (original_number_closed_completions == welltest_state_temp.num_closed_completions()) ) {
453 testWell = false; // this terminates the while loop
454 }
455 }
456
457 // update wellTestState if the well test succeeds
458 if (!welltest_state_temp.well_is_closed(this->name())) {
459 well_test_state.open_well(this->name());
460
461 std::string msg = std::string("well ") + this->name() + std::string(" is re-opened");
462 deferred_logger.info(msg);
463
464 // also reopen completions
465 for (auto& completion : this->well_ecl_.getCompletions()) {
466 if (!welltest_state_temp.completion_is_closed(this->name(), completion.first))
467 well_test_state.open_completion(this->name(), completion.first);
468 }
469 // set the status of the well_state to open
470 ws.open();
471 well_state = well_state_copy;
472 }
473 }
474
475
476
477
478 template<typename TypeTag>
479 bool
481 iterateWellEquations(const Simulator& simulator,
482 const double dt,
483 WellState<Scalar>& well_state,
484 const GroupState<Scalar>& group_state,
485 DeferredLogger& deferred_logger)
486 {
487 const auto& summary_state = simulator.vanguard().summaryState();
488 const auto inj_controls = this->well_ecl_.isInjector() ? this->well_ecl_.injectionControls(summary_state) : Well::InjectionControls(0);
489 const auto prod_controls = this->well_ecl_.isProducer() ? this->well_ecl_.productionControls(summary_state) : Well::ProductionControls(0);
490 bool converged = false;
491 try {
492 // TODO: the following two functions will be refactored to be one to reduce the code duplication
493 if (!this->param_.local_well_solver_control_switching_){
494 converged = this->iterateWellEqWithControl(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
495 } else {
496 if (this->param_.use_implicit_ipr_ && this->well_ecl_.isProducer() && this->wellHasTHPConstraints(summary_state) && (this->well_ecl_.getStatus() == WellStatus::OPEN)) {
497 converged = solveWellWithTHPConstraint(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
498 } else {
499 converged = this->iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
500 }
501 }
502
503 } catch (NumericalProblem& e ) {
504 const std::string msg = "Inner well iterations failed for well " + this->name() + " Treat the well as unconverged. ";
505 deferred_logger.warning("INNER_ITERATION_FAILED", msg);
506 converged = false;
507 }
508 return converged;
509 }
510
511 template<typename TypeTag>
512 bool
515 const double dt,
516 const Well::InjectionControls& inj_controls,
517 const Well::ProductionControls& prod_controls,
518 WellState<Scalar>& well_state,
519 const GroupState<Scalar>& group_state,
520 DeferredLogger& deferred_logger)
521 {
522 const auto& summary_state = simulator.vanguard().summaryState();
523 bool is_operable = true;
524 bool converged = true;
525 auto& ws = well_state.well(this->index_of_well_);
526 // if well is stopped, check if we can reopen
527 if (this->wellIsStopped()) {
528 this->openWell();
529 auto bhp_target = estimateOperableBhp(simulator, dt, well_state, summary_state, deferred_logger);
530 if (!bhp_target.has_value()) {
531 // no intersection with ipr
532 const auto msg = fmt::format("estimateOperableBhp: Did not find operable BHP for well {}", this->name());
533 deferred_logger.debug(msg);
534 is_operable = false;
535 // solve with zero rates
536 solveWellWithZeroRate(simulator, dt, well_state, deferred_logger);
537 this->stopWell();
538 } else {
539 // solve well with the estimated target bhp (or limit)
540 ws.thp = this->getTHPConstraint(summary_state);
541 const Scalar bhp = std::max(bhp_target.value(),
542 static_cast<Scalar>(prod_controls.bhp_limit));
543 solveWellWithBhp(simulator, dt, bhp, well_state, deferred_logger);
544 }
545 }
546 // solve well-equation
547 if (is_operable) {
548 converged = this->iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
549 }
550
551 const bool isThp = ws.production_cmode == Well::ProducerCMode::THP;
552 // check stability of solution under thp-control
553 if (converged && !stoppedOrZeroRateTarget(simulator, well_state, deferred_logger) && isThp) {
554 auto rates = well_state.well(this->index_of_well_).surface_rates;
555 this->adaptRatesForVFP(rates);
556 this->updateIPRImplicit(simulator, well_state, deferred_logger);
557 bool is_stable = WellBhpThpCalculator(*this).isStableSolution(well_state, this->well_ecl_, rates, summary_state);
558 if (!is_stable) {
559 // solution converged to an unstable point!
560 this->operability_status_.use_vfpexplicit = true;
561 auto bhp_stable = WellBhpThpCalculator(*this).estimateStableBhp(well_state, this->well_ecl_, rates, this->getRefDensity(), summary_state);
562 // if we find an intersection with a sufficiently lower bhp, re-solve equations
563 const Scalar reltol = 1e-3;
564 const Scalar cur_bhp = ws.bhp;
565 if (bhp_stable.has_value() && cur_bhp - bhp_stable.value() > cur_bhp*reltol){
566 const auto msg = fmt::format("Well {} converged to an unstable solution, re-solving", this->name());
567 deferred_logger.debug(msg);
568 solveWellWithBhp(simulator, dt, bhp_stable.value(), well_state, deferred_logger);
569 // re-solve with hopefully good initial guess
570 ws.thp = this->getTHPConstraint(summary_state);
571 converged = this->iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
572 }
573 }
574 }
575
576 if (!converged) {
577 // Well did not converge, switch to explicit fractions
578 this->operability_status_.use_vfpexplicit = true;
579 this->openWell();
580 auto bhp_target = estimateOperableBhp(simulator, dt, well_state, summary_state, deferred_logger);
581 if (!bhp_target.has_value()) {
582 // well can't operate using explicit fractions
583 is_operable = false;
584 // solve with zero rate
585 converged = solveWellWithZeroRate(simulator, dt, well_state, deferred_logger);
586 this->stopWell();
587 } else {
588 // solve well with the estimated target bhp (or limit)
589 const Scalar bhp = std::max(bhp_target.value(),
590 static_cast<Scalar>(prod_controls.bhp_limit));
591 solveWellWithBhp(simulator, dt, bhp, well_state, deferred_logger);
592 ws.thp = this->getTHPConstraint(summary_state);
593 converged = this->iterateWellEqWithSwitching(simulator, dt,
594 inj_controls,
595 prod_controls,
596 well_state,
597 group_state,
598 deferred_logger);
599 }
600 }
601 // update operability
602 is_operable = is_operable && !this->wellIsStopped();
603 this->operability_status_.can_obtain_bhp_with_thp_limit = is_operable;
604 this->operability_status_.obey_thp_limit_under_bhp_limit = is_operable;
605 return converged;
606 }
607
608 template<typename TypeTag>
609 std::optional<typename WellInterface<TypeTag>::Scalar>
611 estimateOperableBhp(const Simulator& simulator,
612 const double dt,
613 WellState<Scalar>& well_state,
614 const SummaryState& summary_state,
615 DeferredLogger& deferred_logger)
616 {
617 // Given an unconverged well or closed well, estimate an operable bhp (if any)
618 // Get minimal bhp from vfp-curve
619 Scalar bhp_min = WellBhpThpCalculator(*this).calculateMinimumBhpFromThp(well_state, this->well_ecl_, summary_state, this->getRefDensity());
620 // Solve
621 const bool converged = solveWellWithBhp(simulator, dt, bhp_min, well_state, deferred_logger);
622 if (!converged || this->wellIsStopped()) {
623 return std::nullopt;
624 }
625 this->updateIPRImplicit(simulator, well_state, deferred_logger);
626 auto rates = well_state.well(this->index_of_well_).surface_rates;
627 this->adaptRatesForVFP(rates);
628 return WellBhpThpCalculator(*this).estimateStableBhp(well_state, this->well_ecl_, rates, this->getRefDensity(), summary_state);
629 }
630
631 template<typename TypeTag>
632 bool
634 solveWellWithBhp(const Simulator& simulator,
635 const double dt,
636 const Scalar bhp,
637 WellState<Scalar>& well_state,
638 DeferredLogger& deferred_logger)
639 {
640 // Solve a well using single bhp-constraint (but close if not operable under this)
641 auto group_state = GroupState<Scalar>(); // empty group
642 auto inj_controls = Well::InjectionControls(0);
643 auto prod_controls = Well::ProductionControls(0);
644 auto& ws = well_state.well(this->index_of_well_);
645 auto cmode_inj = ws.injection_cmode;
646 auto cmode_prod = ws.production_cmode;
647 if (this->isInjector()) {
648 inj_controls.addControl(Well::InjectorCMode::BHP);
649 inj_controls.bhp_limit = bhp;
650 inj_controls.cmode = Well::InjectorCMode::BHP;
651 ws.injection_cmode = Well::InjectorCMode::BHP;
652 } else {
653 prod_controls.addControl(Well::ProducerCMode::BHP);
654 prod_controls.bhp_limit = bhp;
655 prod_controls.cmode = Well::ProducerCMode::BHP;
656 ws.production_cmode = Well::ProducerCMode::BHP;
657 }
658 // update well-state
659 ws.bhp = bhp;
660 // solve
661 const bool converged = this->iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger, /*fixed_control*/true);
662 ws.injection_cmode = cmode_inj;
663 ws.production_cmode = cmode_prod;
664 return converged;
665 }
666
667 template<typename TypeTag>
668 bool
670 solveWellWithZeroRate(const Simulator& simulator,
671 const double dt,
672 WellState<Scalar>& well_state,
673 DeferredLogger& deferred_logger)
674 {
675 // Solve a well as stopped
676 const auto well_status_orig = this->wellStatus_;
677 this->stopWell();
678
679 auto group_state = GroupState<Scalar>(); // empty group
680 auto inj_controls = Well::InjectionControls(0);
681 auto prod_controls = Well::ProductionControls(0);
682 const bool converged = this->iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger, /*fixed_control*/true, /*fixed_status*/ true);
683 this->wellStatus_ = well_status_orig;
684 return converged;
685 }
686
687 template<typename TypeTag>
688 bool
690 solveWellForTesting(const Simulator& simulator,
691 WellState<Scalar>& well_state,
692 const GroupState<Scalar>& group_state,
693 DeferredLogger& deferred_logger)
694 {
695 // keep a copy of the original well state
696 const WellState<Scalar> well_state0 = well_state;
697 const double dt = simulator.timeStepSize();
698 const auto& summary_state = simulator.vanguard().summaryState();
699 const bool has_thp_limit = this->wellHasTHPConstraints(summary_state);
700 bool converged;
701 if (has_thp_limit) {
702 well_state.well(this->indexOfWell()).production_cmode = Well::ProducerCMode::THP;
703 converged = gliftBeginTimeStepWellTestIterateWellEquations(
704 simulator, dt, well_state, group_state, deferred_logger);
705 }
706 else {
707 well_state.well(this->indexOfWell()).production_cmode = Well::ProducerCMode::BHP;
708 converged = iterateWellEquations(simulator, dt, well_state, group_state, deferred_logger);
709 }
710 if (converged) {
711 deferred_logger.debug("WellTest: Well equation for well " + this->name() + " converged");
712 return true;
713 }
714 const int max_iter = param_.max_welleq_iter_;
715 deferred_logger.debug("WellTest: Well equation for well " + this->name() + " failed converging in "
716 + std::to_string(max_iter) + " iterations");
717 well_state = well_state0;
718 return false;
719 }
720
721
722 template<typename TypeTag>
723 void
725 solveWellEquation(const Simulator& simulator,
726 WellState<Scalar>& well_state,
727 const GroupState<Scalar>& group_state,
728 DeferredLogger& deferred_logger)
729 {
730 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
731 return;
732
733 // keep a copy of the original well state
734 const WellState<Scalar> well_state0 = well_state;
735 const double dt = simulator.timeStepSize();
736 bool converged = iterateWellEquations(simulator, dt, well_state, group_state, deferred_logger);
737
738 // Newly opened wells with THP control sometimes struggles to
739 // converge due to bad initial guess. Or due to the simple fact
740 // that the well needs to change to another control.
741 // We therefore try to solve the well with BHP control to get
742 // an better initial guess.
743 // If the well is supposed to operate under THP control
744 // "updateWellControl" will switch it back to THP later.
745 if (!converged) {
746 auto& ws = well_state.well(this->indexOfWell());
747 bool thp_control = false;
748 if (this->well_ecl_.isInjector()) {
749 thp_control = ws.injection_cmode == Well::InjectorCMode::THP;
750 if (thp_control) {
751 ws.injection_cmode = Well::InjectorCMode::BHP;
752 this->well_control_log_.push_back(WellInjectorCMode2String(Well::InjectorCMode::THP));
753 }
754 } else {
755 thp_control = ws.production_cmode == Well::ProducerCMode::THP;
756 if (thp_control) {
757 ws.production_cmode = Well::ProducerCMode::BHP;
758 this->well_control_log_.push_back(WellProducerCMode2String(Well::ProducerCMode::THP));
759 }
760 }
761 if (thp_control) {
762 const std::string msg = std::string("The newly opened well ") + this->name()
763 + std::string(" with THP control did not converge during inner iterations, we try again with bhp control");
764 deferred_logger.debug(msg);
765 converged = this->iterateWellEquations(simulator, dt, well_state, group_state, deferred_logger);
766 }
767 }
768
769 if (!converged) {
770 const int max_iter = param_.max_welleq_iter_;
771 deferred_logger.debug("Compute initial well solution for well " + this->name() + ". Failed to converge in "
772 + std::to_string(max_iter) + " iterations");
773 well_state = well_state0;
774 }
775 }
776
777
778
779 template <typename TypeTag>
780 void
782 assembleWellEq(const Simulator& simulator,
783 const double dt,
784 WellState<Scalar>& well_state,
785 const GroupState<Scalar>& group_state,
786 DeferredLogger& deferred_logger)
787 {
788 prepareWellBeforeAssembling(simulator, dt, well_state, group_state, deferred_logger);
789 assembleWellEqWithoutIteration(simulator, dt, well_state, group_state, deferred_logger);
790 }
791
792
793
794 template <typename TypeTag>
795 void
798 const double dt,
799 WellState<Scalar>& well_state,
800 const GroupState<Scalar>& group_state,
801 DeferredLogger& deferred_logger)
802 {
803 const auto& summary_state = simulator.vanguard().summaryState();
804 const auto inj_controls = this->well_ecl_.isInjector() ? this->well_ecl_.injectionControls(summary_state) : Well::InjectionControls(0);
805 const auto prod_controls = this->well_ecl_.isProducer() ? this->well_ecl_.productionControls(summary_state) : Well::ProductionControls(0);
806 // TODO: the reason to have inj_controls and prod_controls in the arguments, is that we want to change the control used for the well functions
807 // TODO: maybe we can use std::optional or pointers to simplify here
808 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
809 }
810
811
812
813 template<typename TypeTag>
814 void
817 const double dt,
818 WellState<Scalar>& well_state,
819 const GroupState<Scalar>& group_state,
820 DeferredLogger& deferred_logger)
821 {
822 const bool old_well_operable = this->operability_status_.isOperableAndSolvable();
823
824 if (param_.check_well_operability_iter_)
825 checkWellOperability(simulator, well_state, deferred_logger);
826
827 // only use inner well iterations for the first newton iterations.
828 const int iteration_idx = simulator.model().newtonMethod().numIterations();
829 if (iteration_idx < param_.max_niter_inner_well_iter_ || this->well_ecl_.isMultiSegment()) {
830 this->operability_status_.solvable = true;
831 bool converged = this->iterateWellEquations(simulator, dt, well_state, group_state, deferred_logger);
832
833 // unsolvable wells are treated as not operable and will not be solved for in this iteration.
834 if (!converged) {
835 if (param_.shut_unsolvable_wells_)
836 this->operability_status_.solvable = false;
837 }
838 }
839 if (this->operability_status_.has_negative_potentials) {
840 auto well_state_copy = well_state;
841 std::vector<Scalar> potentials;
842 try {
843 computeWellPotentials(simulator, well_state_copy, potentials, deferred_logger);
844 } catch (const std::exception& e) {
845 const std::string msg = fmt::format("well {}: computeWellPotentials() failed "
846 "during attempt to recompute potentials for well: ",
847 this->name(), e.what());
848 deferred_logger.info(msg);
849 this->operability_status_.has_negative_potentials = true;
850 }
851 auto& ws = well_state.well(this->indexOfWell());
852 const int np = well_state.numPhases();
853 for (int p = 0; p < np; ++p) {
854 ws.well_potentials[p] = std::max(Scalar{0.0}, potentials[p]);
855 }
856 }
857 this->changed_to_open_this_step_ = false;
858 const bool well_operable = this->operability_status_.isOperableAndSolvable();
859
860 if (!well_operable && old_well_operable) {
861 deferred_logger.info(" well " + this->name() + " gets STOPPED during iteration ");
862 this->stopWell();
863 changed_to_stopped_this_step_ = true;
864 } else if (well_operable && !old_well_operable) {
865 deferred_logger.info(" well " + this->name() + " gets REVIVED during iteration ");
866 this->openWell();
867 changed_to_stopped_this_step_ = false;
868 this->changed_to_open_this_step_ = true;
869 }
870 }
871
872 template<typename TypeTag>
873 void
875 {
876 if(!this->isOperableAndSolvable() && !this->wellIsStopped())
877 return;
878
879 for (int perfIdx = 0; perfIdx < this->number_of_perforations_; ++perfIdx) {
880 if (this->cells()[perfIdx] == cellIdx) {
881 for (int i = 0; i < RateVector::dimension; ++i) {
882 rates[i] += connectionRates_[perfIdx][i];
883 }
884 }
885 }
886 }
887
888 template<typename TypeTag>
891 {
892 for (int perfIdx = 0; perfIdx < this->number_of_perforations_; ++perfIdx) {
893 if (this->cells()[perfIdx] == cellIdx) {
894 const unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
895 return connectionRates_[perfIdx][activeCompIdx].value();
896 }
897 }
898 // this is not thread safe
899 OPM_THROW(std::invalid_argument, "The well with name " + this->name()
900 + " does not perforate cell " + std::to_string(cellIdx));
901 return 0.0;
902 }
903
904
905
906
907 template<typename TypeTag>
908 void
910 checkWellOperability(const Simulator& simulator,
911 const WellState<Scalar>& well_state,
912 DeferredLogger& deferred_logger)
913 {
914
915 if (!param_.check_well_operability_) {
916 return;
917 }
918
919 if (this->wellIsStopped() && !changed_to_stopped_this_step_) {
920 return;
921 }
922
923 updateWellOperability(simulator, well_state, deferred_logger);
924 if (!this->operability_status_.isOperableAndSolvable()) {
925 this->operability_status_.use_vfpexplicit = true;
926 deferred_logger.debug("EXPLICIT_LOOKUP_VFP",
927 "well not operable, trying with explicit vfp lookup: " + this->name());
928 updateWellOperability(simulator, well_state, deferred_logger);
929 }
930 }
931
932 template<typename TypeTag>
933 bool
936 const double dt,
937 WellState<Scalar>& well_state,
938 const GroupState<Scalar>& group_state,
939 DeferredLogger& deferred_logger)
940 {
941 const auto& well_name = this->name();
942 assert(this->wellHasTHPConstraints(simulator.vanguard().summaryState()));
943 const auto& schedule = simulator.vanguard().schedule();
944 auto report_step_idx = simulator.episodeIndex();
945 const auto& glo = schedule.glo(report_step_idx);
946 if(glo.active() && glo.has_well(well_name)) {
947 const auto increment = glo.gaslift_increment();
948 auto alq = well_state.getALQ(well_name);
949 bool converged;
950 while (alq > 0) {
951 well_state.setALQ(well_name, alq);
952 if ((converged =
953 iterateWellEquations(simulator, dt, well_state, group_state, deferred_logger)))
954 {
955 return converged;
956 }
957 alq -= increment;
958 }
959 return false;
960 }
961 else {
962 return iterateWellEquations(simulator, dt, well_state, group_state, deferred_logger);
963 }
964 }
965
966 template<typename TypeTag>
967 void
970 WellState<Scalar>& well_state,
971 DeferredLogger& deferred_logger)
972 {
973 const auto& summary_state = simulator.vanguard().summaryState();
974 const auto& well_name = this->name();
975 if (!this->wellHasTHPConstraints(summary_state)) {
976 const std::string msg = fmt::format("GLIFT WTEST: Well {} does not have THP constraints", well_name);
977 deferred_logger.info(msg);
978 return;
979 }
980 const auto& schedule = simulator.vanguard().schedule();
981 const auto report_step_idx = simulator.episodeIndex();
982 const auto& glo = schedule.glo(report_step_idx);
983 if (!glo.has_well(well_name)) {
984 const std::string msg = fmt::format(
985 "GLIFT WTEST: Well {} : Gas Lift not activated: "
986 "WLIFTOPT is probably missing. Skipping.", well_name);
987 deferred_logger.info(msg);
988 return;
989 }
990 const auto& gl_well = glo.well(well_name);
991 auto& max_alq_optional = gl_well.max_rate();
992 Scalar max_alq;
993 if (max_alq_optional) {
994 max_alq = *max_alq_optional;
995 }
996 else {
997 const auto& well_ecl = this->wellEcl();
998 const auto& controls = well_ecl.productionControls(summary_state);
999 const auto& table = this->vfpProperties()->getProd()->getTable(controls.vfp_table_number);
1000 const auto& alq_values = table.getALQAxis();
1001 max_alq = alq_values.back();
1002 }
1003 well_state.setALQ(well_name, max_alq);
1004 const std::string msg = fmt::format(
1005 "GLIFT WTEST: Well {} : Setting ALQ to max value: {}",
1006 well_name, max_alq);
1007 deferred_logger.info(msg);
1008 }
1009
1010 template<typename TypeTag>
1011 void
1013 updateWellOperability(const Simulator& simulator,
1014 const WellState<Scalar>& well_state,
1015 DeferredLogger& deferred_logger)
1016 {
1017 if (this->param_.local_well_solver_control_switching_) {
1018 const bool success = updateWellOperabilityFromWellEq(simulator, well_state, deferred_logger);
1019 if (success) {
1020 return;
1021 } else {
1022 deferred_logger.debug("Operability check using well equations did not converge for well "
1023 + this->name() + ", reverting to classical approach." );
1024 }
1025 }
1026 this->operability_status_.resetOperability();
1027
1028 bool thp_controlled = this->isInjector() ? well_state.well(this->index_of_well_).injection_cmode == Well::InjectorCMode::THP:
1029 well_state.well(this->index_of_well_).production_cmode == Well::ProducerCMode::THP;
1030 bool bhp_controlled = this->isInjector() ? well_state.well(this->index_of_well_).injection_cmode == Well::InjectorCMode::BHP:
1031 well_state.well(this->index_of_well_).production_cmode == Well::ProducerCMode::BHP;
1032
1033 // Operability checking is not free
1034 // Only check wells under BHP and THP control
1035 bool check_thp = thp_controlled || this->operability_status_.thp_limit_violated_but_not_switched;
1036 if (check_thp || bhp_controlled) {
1037 updateIPR(simulator, deferred_logger);
1038 checkOperabilityUnderBHPLimit(well_state, simulator, deferred_logger);
1039 }
1040 // we do some extra checking for wells under THP control.
1041 if (check_thp) {
1042 checkOperabilityUnderTHPLimit(simulator, well_state, deferred_logger);
1043 }
1044 }
1045
1046 template<typename TypeTag>
1047 bool
1050 const WellState<Scalar>& well_state,
1051 DeferredLogger& deferred_logger)
1052 {
1053 // only makes sense if we're using this parameter is true
1054 assert(this->param_.local_well_solver_control_switching_);
1055 this->operability_status_.resetOperability();
1056 WellState<Scalar> well_state_copy = well_state;
1057 const auto& group_state = simulator.problem().wellModel().groupState();
1058 const double dt = simulator.timeStepSize();
1059 // equations should be converged at this stage, so only one it is needed
1060 bool converged = iterateWellEquations(simulator, dt, well_state_copy, group_state, deferred_logger);
1061 return converged;
1062 }
1063
1064 template<typename TypeTag>
1065 void
1067 updateWellStateWithTarget(const Simulator& simulator,
1068 const GroupState<Scalar>& group_state,
1069 WellState<Scalar>& well_state,
1070 DeferredLogger& deferred_logger) const
1071 {
1072
1073 // only bhp and wellRates are used to initilize the primaryvariables for standard wells
1074 const auto& well = this->well_ecl_;
1075 const int well_index = this->index_of_well_;
1076 auto& ws = well_state.well(well_index);
1077 const auto& pu = this->phaseUsage();
1078 const int np = well_state.numPhases();
1079 const auto& summaryState = simulator.vanguard().summaryState();
1080 const auto& schedule = simulator.vanguard().schedule();
1081
1082 if (this->wellIsStopped()) {
1083 for (int p = 0; p<np; ++p) {
1084 ws.surface_rates[p] = 0;
1085 }
1086 ws.thp = 0;
1087 return;
1088 }
1089
1090 if (this->isInjector() )
1091 {
1092 const auto& controls = well.injectionControls(summaryState);
1093
1094 InjectorType injectorType = controls.injector_type;
1095 int phasePos;
1096 switch (injectorType) {
1097 case InjectorType::WATER:
1098 {
1099 phasePos = pu.phase_pos[BlackoilPhases::Aqua];
1100 break;
1101 }
1102 case InjectorType::OIL:
1103 {
1104 phasePos = pu.phase_pos[BlackoilPhases::Liquid];
1105 break;
1106 }
1107 case InjectorType::GAS:
1108 {
1109 phasePos = pu.phase_pos[BlackoilPhases::Vapour];
1110 break;
1111 }
1112 default:
1113 OPM_DEFLOG_THROW(std::runtime_error, "Expected WATER, OIL or GAS as type for injectors " + this->name(), deferred_logger );
1114 }
1115
1116 const auto current = ws.injection_cmode;
1117
1118 switch (current) {
1119 case Well::InjectorCMode::RATE:
1120 {
1121 ws.surface_rates[phasePos] = (1.0 - this->rsRvInj()) * controls.surface_rate;
1122 if(this->rsRvInj() > 0) {
1123 if (injectorType == InjectorType::OIL && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1124 ws.surface_rates[pu.phase_pos[BlackoilPhases::Vapour]] = controls.surface_rate * this->rsRvInj();
1125 } else if (injectorType == InjectorType::GAS && FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1126 ws.surface_rates[pu.phase_pos[BlackoilPhases::Liquid]] = controls.surface_rate * this->rsRvInj();
1127 } else {
1128 OPM_DEFLOG_THROW(std::runtime_error, "Expected OIL or GAS as type for injectors when RS/RV (item 10) is non-zero " + this->name(), deferred_logger );
1129 }
1130 }
1131 break;
1132 }
1133
1134 case Well::InjectorCMode::RESV:
1135 {
1136 std::vector<Scalar> convert_coeff(this->number_of_phases_, 1.0);
1137 this->rateConverter_.calcCoeff(/*fipreg*/ 0, this->pvtRegionIdx_, convert_coeff);
1138 const Scalar coeff = convert_coeff[phasePos];
1139 ws.surface_rates[phasePos] = controls.reservoir_rate/coeff;
1140 break;
1141 }
1142
1143 case Well::InjectorCMode::THP:
1144 {
1145 auto rates = ws.surface_rates;
1146 Scalar bhp = WellBhpThpCalculator(*this).calculateBhpFromThp(well_state,
1147 rates,
1148 well,
1149 summaryState,
1150 this->getRefDensity(),
1151 deferred_logger);
1152 ws.bhp = bhp;
1153 ws.thp = this->getTHPConstraint(summaryState);
1154
1155 // if the total rates are negative or zero
1156 // we try to provide a better intial well rate
1157 // using the well potentials
1158 Scalar total_rate = std::accumulate(rates.begin(), rates.end(), 0.0);
1159 if (total_rate <= 0.0)
1160 ws.surface_rates = ws.well_potentials;
1161
1162 break;
1163 }
1164 case Well::InjectorCMode::BHP:
1165 {
1166 ws.bhp = controls.bhp_limit;
1167 Scalar total_rate = 0.0;
1168 for (int p = 0; p<np; ++p) {
1169 total_rate += ws.surface_rates[p];
1170 }
1171 // if the total rates are negative or zero
1172 // we try to provide a better intial well rate
1173 // using the well potentials
1174 if (total_rate <= 0.0)
1175 ws.surface_rates = ws.well_potentials;
1176
1177 break;
1178 }
1179 case Well::InjectorCMode::GRUP:
1180 {
1181 assert(well.isAvailableForGroupControl());
1182 const auto& group = schedule.getGroup(well.groupName(), this->currentStep());
1183 const Scalar efficiencyFactor = well.getEfficiencyFactor();
1184 std::optional<Scalar> target =
1185 this->getGroupInjectionTargetRate(group,
1186 well_state,
1187 group_state,
1188 schedule,
1189 summaryState,
1190 injectorType,
1191 efficiencyFactor,
1192 deferred_logger);
1193 if (target)
1194 ws.surface_rates[phasePos] = *target;
1195 break;
1196 }
1197 case Well::InjectorCMode::CMODE_UNDEFINED:
1198 {
1199 OPM_DEFLOG_THROW(std::runtime_error, "Well control must be specified for well " + this->name(), deferred_logger );
1200 }
1201
1202 }
1203 // for wells with zero injection rate, if we assign exactly zero rate,
1204 // we will have to assume some trivial composition in the wellbore.
1205 // here, we use some small value (about 0.01 m^3/day ~= 1.e-7) to initialize
1206 // the zero rate target, then we can use to retain the composition information
1207 // within the wellbore from the previous result, and hopefully it is a good
1208 // initial guess for the zero rate target.
1209 ws.surface_rates[phasePos] = std::max(Scalar{1.e-7}, ws.surface_rates[phasePos]);
1210
1211 if (ws.bhp == 0.) {
1212 ws.bhp = controls.bhp_limit;
1213 }
1214 }
1215 //Producer
1216 else
1217 {
1218 const auto current = ws.production_cmode;
1219 const auto& controls = well.productionControls(summaryState);
1220 switch (current) {
1221 case Well::ProducerCMode::ORAT:
1222 {
1223 Scalar current_rate = -ws.surface_rates[ pu.phase_pos[Oil] ];
1224 // for trivial rates or opposite direction we don't just scale the rates
1225 // but use either the potentials or the mobility ratio to initial the well rates
1226 if (current_rate > 0.0) {
1227 for (int p = 0; p<np; ++p) {
1228 ws.surface_rates[p] *= controls.oil_rate/current_rate;
1229 }
1230 } else {
1231 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1232 double control_fraction = fractions[pu.phase_pos[Oil]];
1233 if (control_fraction != 0.0) {
1234 for (int p = 0; p<np; ++p) {
1235 ws.surface_rates[p] = - fractions[p] * controls.oil_rate/control_fraction;
1236 }
1237 }
1238 }
1239 break;
1240 }
1241 case Well::ProducerCMode::WRAT:
1242 {
1243 Scalar current_rate = -ws.surface_rates[ pu.phase_pos[Water] ];
1244 // for trivial rates or opposite direction we don't just scale the rates
1245 // but use either the potentials or the mobility ratio to initial the well rates
1246 if (current_rate > 0.0) {
1247 for (int p = 0; p<np; ++p) {
1248 ws.surface_rates[p] *= controls.water_rate/current_rate;
1249 }
1250 } else {
1251 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1252 const Scalar control_fraction = fractions[pu.phase_pos[Water]];
1253 if (control_fraction != 0.0) {
1254 for (int p = 0; p<np; ++p) {
1255 ws.surface_rates[p] = - fractions[p] * controls.water_rate / control_fraction;
1256 }
1257 }
1258 }
1259 break;
1260 }
1261 case Well::ProducerCMode::GRAT:
1262 {
1263 Scalar current_rate = -ws.surface_rates[pu.phase_pos[Gas] ];
1264 // or trivial rates or opposite direction we don't just scale the rates
1265 // but use either the potentials or the mobility ratio to initial the well rates
1266 if (current_rate > 0.0) {
1267 for (int p = 0; p<np; ++p) {
1268 ws.surface_rates[p] *= controls.gas_rate/current_rate;
1269 }
1270 } else {
1271 const std::vector<Scalar > fractions = initialWellRateFractions(simulator, well_state);
1272 const Scalar control_fraction = fractions[pu.phase_pos[Gas]];
1273 if (control_fraction != 0.0) {
1274 for (int p = 0; p<np; ++p) {
1275 ws.surface_rates[p] = - fractions[p] * controls.gas_rate / control_fraction;
1276 }
1277 }
1278 }
1279
1280 break;
1281
1282 }
1283 case Well::ProducerCMode::LRAT:
1284 {
1285 Scalar current_rate = - ws.surface_rates[ pu.phase_pos[Water] ]
1286 - ws.surface_rates[ pu.phase_pos[Oil] ];
1287 // or trivial rates or opposite direction we don't just scale the rates
1288 // but use either the potentials or the mobility ratio to initial the well rates
1289 if (current_rate > 0.0) {
1290 for (int p = 0; p<np; ++p) {
1291 ws.surface_rates[p] *= controls.liquid_rate/current_rate;
1292 }
1293 } else {
1294 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1295 const Scalar control_fraction = fractions[pu.phase_pos[Water]] + fractions[pu.phase_pos[Oil]];
1296 if (control_fraction != 0.0) {
1297 for (int p = 0; p<np; ++p) {
1298 ws.surface_rates[p] = - fractions[p] * controls.liquid_rate / control_fraction;
1299 }
1300 }
1301 }
1302 break;
1303 }
1304 case Well::ProducerCMode::CRAT:
1305 {
1306 OPM_DEFLOG_THROW(std::runtime_error,
1307 fmt::format("CRAT control not supported, well {}", this->name()),
1308 deferred_logger);
1309 }
1310 case Well::ProducerCMode::RESV:
1311 {
1312 std::vector<Scalar> convert_coeff(this->number_of_phases_, 1.0);
1313 this->rateConverter_.calcCoeff(/*fipreg*/ 0, this->pvtRegionIdx_, ws.surface_rates, convert_coeff);
1314 Scalar total_res_rate = 0.0;
1315 for (int p = 0; p<np; ++p) {
1316 total_res_rate -= ws.surface_rates[p] * convert_coeff[p];
1317 }
1318 if (controls.prediction_mode) {
1319 // or trivial rates or opposite direction we don't just scale the rates
1320 // but use either the potentials or the mobility ratio to initial the well rates
1321 if (total_res_rate > 0.0) {
1322 for (int p = 0; p<np; ++p) {
1323 ws.surface_rates[p] *= controls.resv_rate/total_res_rate;
1324 }
1325 } else {
1326 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1327 for (int p = 0; p<np; ++p) {
1328 ws.surface_rates[p] = - fractions[p] * controls.resv_rate / convert_coeff[p];
1329 }
1330 }
1331 } else {
1332 std::vector<Scalar> hrates(this->number_of_phases_,0.);
1333 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
1334 hrates[pu.phase_pos[Water]] = controls.water_rate;
1335 }
1336 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1337 hrates[pu.phase_pos[Oil]] = controls.oil_rate;
1338 }
1339 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1340 hrates[pu.phase_pos[Gas]] = controls.gas_rate;
1341 }
1342 std::vector<Scalar> hrates_resv(this->number_of_phases_,0.);
1343 this->rateConverter_.calcReservoirVoidageRates(/*fipreg*/ 0, this->pvtRegionIdx_, hrates, hrates_resv);
1344 Scalar target = std::accumulate(hrates_resv.begin(), hrates_resv.end(), 0.0);
1345 // or trivial rates or opposite direction we don't just scale the rates
1346 // but use either the potentials or the mobility ratio to initial the well rates
1347 if (total_res_rate > 0.0) {
1348 for (int p = 0; p<np; ++p) {
1349 ws.surface_rates[p] *= target/total_res_rate;
1350 }
1351 } else {
1352 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1353 for (int p = 0; p<np; ++p) {
1354 ws.surface_rates[p] = - fractions[p] * target / convert_coeff[p];
1355 }
1356 }
1357 }
1358 break;
1359 }
1360 case Well::ProducerCMode::BHP:
1361 {
1362 ws.bhp = controls.bhp_limit;
1363 Scalar total_rate = 0.0;
1364 for (int p = 0; p<np; ++p) {
1365 total_rate -= ws.surface_rates[p];
1366 }
1367 // if the total rates are negative or zero
1368 // we try to provide a better intial well rate
1369 // using the well potentials
1370 if (total_rate <= 0.0){
1371 for (int p = 0; p<np; ++p) {
1372 ws.surface_rates[p] = -ws.well_potentials[p];
1373 }
1374 }
1375 break;
1376 }
1377 case Well::ProducerCMode::THP:
1378 {
1379 const bool update_success = updateWellStateWithTHPTargetProd(simulator, well_state, deferred_logger);
1380
1381 if (!update_success) {
1382 // the following is the original way of initializing well state with THP constraint
1383 // keeping it for robust reason in case that it fails to get a bhp value with THP constraint
1384 // more sophisticated design might be needed in the future
1385 auto rates = ws.surface_rates;
1386 this->adaptRatesForVFP(rates);
1388 well_state, rates, well, summaryState, this->getRefDensity(), deferred_logger);
1389 ws.bhp = bhp;
1390 ws.thp = this->getTHPConstraint(summaryState);
1391 // if the total rates are negative or zero
1392 // we try to provide a better initial well rate
1393 // using the well potentials
1394 const Scalar total_rate = -std::accumulate(rates.begin(), rates.end(), 0.0);
1395 if (total_rate <= 0.0) {
1396 for (int p = 0; p < this->number_of_phases_; ++p) {
1397 ws.surface_rates[p] = -ws.well_potentials[p];
1398 }
1399 }
1400 }
1401 break;
1402 }
1403 case Well::ProducerCMode::GRUP:
1404 {
1405 assert(well.isAvailableForGroupControl());
1406 const auto& group = schedule.getGroup(well.groupName(), this->currentStep());
1407 const Scalar efficiencyFactor = well.getEfficiencyFactor();
1408 Scalar scale = this->getGroupProductionTargetRate(group,
1409 well_state,
1410 group_state,
1411 schedule,
1412 summaryState,
1413 efficiencyFactor,
1414 deferred_logger);
1415
1416 // we don't want to scale with zero and get zero rates.
1417 if (scale > 0) {
1418 for (int p = 0; p<np; ++p) {
1419 ws.surface_rates[p] *= scale;
1420 }
1421 ws.trivial_target = false;
1422 } else {
1423 ws.trivial_target = true;
1424 }
1425 break;
1426 }
1427 case Well::ProducerCMode::CMODE_UNDEFINED:
1429 {
1430 OPM_DEFLOG_THROW(std::runtime_error, "Well control must be specified for well " + this->name() , deferred_logger);
1431 break;
1432 }
1433 } // end of switch
1434
1435 if (ws.bhp == 0.) {
1436 ws.bhp = controls.bhp_limit;
1437 }
1438 }
1439 }
1440
1441 template<typename TypeTag>
1442 bool
1444 wellUnderZeroRateTarget(const Simulator& simulator,
1445 const WellState<Scalar>& well_state,
1446 DeferredLogger& deferred_logger) const
1447 {
1448 // Check if well is under zero rate control, either directly or from group
1449 const bool isGroupControlled = this->wellUnderGroupControl(well_state.well(this->index_of_well_));
1450 if (!isGroupControlled) {
1451 // well is not under group control, check "individual" version
1452 const auto& summaryState = simulator.vanguard().summaryState();
1453 return this->wellUnderZeroRateTargetIndividual(summaryState, well_state);
1454 } else {
1455 return this->wellUnderZeroGroupRateTarget(simulator, well_state, deferred_logger, isGroupControlled);
1456 }
1457 }
1458
1459 template <typename TypeTag>
1460 bool
1462 const WellState<Scalar>& well_state,
1463 DeferredLogger& deferred_logger,
1464 const std::optional<bool> group_control) const
1465 {
1466 // Check if well is under zero rate target from group
1467 const bool isGroupControlled = group_control.value_or(this->wellUnderGroupControl(well_state.well(this->index_of_well_)));
1468 if (isGroupControlled) {
1469 const auto& summaryState = simulator.vanguard().summaryState();
1470 const auto& group_state = simulator.problem().wellModel().groupState();
1471 const auto& schedule = simulator.vanguard().schedule();
1472 return this->zeroGroupRateTarget(summaryState, schedule, well_state, group_state, deferred_logger);
1473 }
1474 return false;
1475 }
1476
1477 template<typename TypeTag>
1478 bool
1480 stoppedOrZeroRateTarget(const Simulator& simulator,
1481 const WellState<Scalar>& well_state,
1482 DeferredLogger& deferred_logger) const
1483 {
1484 // Check if well is stopped or under zero rate control, either
1485 // directly or from group.
1486 return this->wellIsStopped()
1487 || this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger);
1488 }
1489
1490 template<typename TypeTag>
1491 std::vector<typename WellInterface<TypeTag>::Scalar>
1493 initialWellRateFractions(const Simulator& simulator,
1494 const WellState<Scalar>& well_state) const
1495 {
1496 const int np = this->number_of_phases_;
1497 std::vector<Scalar> scaling_factor(np);
1498 const auto& ws = well_state.well(this->index_of_well_);
1499
1500 Scalar total_potentials = 0.0;
1501 for (int p = 0; p<np; ++p) {
1502 total_potentials += ws.well_potentials[p];
1503 }
1504 if (total_potentials > 0) {
1505 for (int p = 0; p<np; ++p) {
1506 scaling_factor[p] = ws.well_potentials[p] / total_potentials;
1507 }
1508 return scaling_factor;
1509 }
1510 // if we don't have any potentials we weight it using the mobilites
1511 // We only need approximation so we don't bother with the vapporized oil and dissolved gas
1512 Scalar total_tw = 0;
1513 const int nperf = this->number_of_perforations_;
1514 for (int perf = 0; perf < nperf; ++perf) {
1515 total_tw += this->well_index_[perf];
1516 }
1517 for (int perf = 0; perf < nperf; ++perf) {
1518 const int cell_idx = this->well_cells_[perf];
1519 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
1520 const auto& fs = intQuants.fluidState();
1521 const Scalar well_tw_fraction = this->well_index_[perf] / total_tw;
1522 Scalar total_mobility = 0.0;
1523 for (int p = 0; p < np; ++p) {
1524 int modelPhaseIdx = this->flowPhaseToModelPhaseIdx(p);
1525 total_mobility += fs.invB(modelPhaseIdx).value() * intQuants.mobility(modelPhaseIdx).value();
1526 }
1527 for (int p = 0; p < np; ++p) {
1528 int modelPhaseIdx = this->flowPhaseToModelPhaseIdx(p);
1529 scaling_factor[p] += well_tw_fraction * fs.invB(modelPhaseIdx).value() * intQuants.mobility(modelPhaseIdx).value() / total_mobility;
1530 }
1531 }
1532 return scaling_factor;
1533 }
1534
1535
1536
1537 template <typename TypeTag>
1538 void
1540 updateWellStateRates(const Simulator& simulator,
1541 WellState<Scalar>& well_state,
1542 DeferredLogger& deferred_logger) const
1543 {
1544 // Check if the rates of this well only are single-phase, do nothing
1545 // if more than one nonzero rate.
1546 auto& ws = well_state.well(this->index_of_well_);
1547 int nonzero_rate_index = -1;
1548 const Scalar floating_point_error_epsilon = 1e-14;
1549 for (int p = 0; p < this->number_of_phases_; ++p) {
1550 if (std::abs(ws.surface_rates[p]) > floating_point_error_epsilon) {
1551 if (nonzero_rate_index == -1) {
1552 nonzero_rate_index = p;
1553 } else {
1554 // More than one nonzero rate.
1555 return;
1556 }
1557 }
1558 }
1559
1560 // Calculate the rates that follow from the current primary variables.
1561 std::vector<Scalar> well_q_s = computeCurrentWellRates(simulator, deferred_logger);
1562
1563 if (nonzero_rate_index == -1) {
1564 // No nonzero rates.
1565 // Use the computed rate directly
1566 for (int p = 0; p < this->number_of_phases_; ++p) {
1567 ws.surface_rates[p] = well_q_s[this->flowPhaseToModelCompIdx(p)];
1568 }
1569 return;
1570 }
1571
1572 // Set the currently-zero phase flows to be nonzero in proportion to well_q_s.
1573 const Scalar initial_nonzero_rate = ws.surface_rates[nonzero_rate_index];
1574 const int comp_idx_nz = this->flowPhaseToModelCompIdx(nonzero_rate_index);
1575 if (std::abs(well_q_s[comp_idx_nz]) > floating_point_error_epsilon) {
1576 for (int p = 0; p < this->number_of_phases_; ++p) {
1577 if (p != nonzero_rate_index) {
1578 const int comp_idx = this->flowPhaseToModelCompIdx(p);
1579 Scalar& rate = ws.surface_rates[p];
1580 rate = (initial_nonzero_rate / well_q_s[comp_idx_nz]) * (well_q_s[comp_idx]);
1581 }
1582 }
1583 }
1584 }
1585
1586 template <typename TypeTag>
1587 std::vector<typename WellInterface<TypeTag>::Scalar>
1589 wellIndex(const int perf,
1590 const IntensiveQuantities& intQuants,
1591 const Scalar trans_mult,
1592 const SingleWellState<Scalar>& ws) const
1593 {
1594 // Add a Forchheimer term to the gas phase CTF if the run uses
1595 // either of the WDFAC or the WDFACCOR keywords.
1596
1597 auto wi = std::vector<Scalar>
1598 (this->num_components_, this->well_index_[perf] * trans_mult);
1599
1600 if constexpr (! Indices::gasEnabled) {
1601 return wi;
1602 }
1603
1604 const auto& wdfac = this->well_ecl_.getWDFAC();
1605
1606 if (! wdfac.useDFactor() || (this->well_index_[perf] == 0.0)) {
1607 return wi;
1608 }
1609
1610 const Scalar d = this->computeConnectionDFactor(perf, intQuants, ws);
1611 if (d < 1.0e-15) {
1612 return wi;
1613 }
1614
1615 // Solve quadratic equations for connection rates satisfying the ipr and the flow-dependent skin.
1616 // If more than one solution, pick the one corresponding to lowest absolute rate (smallest skin).
1617 const auto& connection = this->well_ecl_.getConnections()[ws.perf_data.ecl_index[perf]];
1618 const Scalar Kh = connection.Kh();
1619 const Scalar scaling = 3.141592653589 * Kh * connection.wpimult();
1620 const unsigned gas_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
1621
1622 const Scalar connection_pressure = ws.perf_data.pressure[perf];
1623 const Scalar cell_pressure = getValue(intQuants.fluidState().pressure(FluidSystem::gasPhaseIdx));
1624 const Scalar drawdown = cell_pressure - connection_pressure;
1625 const Scalar invB = getValue(intQuants.fluidState().invB(FluidSystem::gasPhaseIdx));
1626 const Scalar mob_g = getValue(intQuants.mobility(FluidSystem::gasPhaseIdx)) * invB;
1627 const Scalar a = d;
1628 const Scalar b = 2*scaling/wi[gas_comp_idx];
1629 const Scalar c = -2*scaling*mob_g*drawdown;
1630
1631 Scalar consistent_Q = -1.0e20;
1632 // Find and check negative solutions (a --> -a)
1633 const Scalar r2n = b*b + 4*a*c;
1634 if (r2n >= 0) {
1635 const Scalar rn = std::sqrt(r2n);
1636 const Scalar xn1 = (b-rn)*0.5/a;
1637 if (xn1 <= 0) {
1638 consistent_Q = xn1;
1639 }
1640 const Scalar xn2 = (b+rn)*0.5/a;
1641 if (xn2 <= 0 && xn2 > consistent_Q) {
1642 consistent_Q = xn2;
1643 }
1644 }
1645 // Find and check positive solutions
1646 consistent_Q *= -1;
1647 const Scalar r2p = b*b - 4*a*c;
1648 if (r2p >= 0) {
1649 const Scalar rp = std::sqrt(r2p);
1650 const Scalar xp1 = (rp-b)*0.5/a;
1651 if (xp1 > 0 && xp1 < consistent_Q) {
1652 consistent_Q = xp1;
1653 }
1654 const Scalar xp2 = -(rp+b)*0.5/a;
1655 if (xp2 > 0 && xp2 < consistent_Q) {
1656 consistent_Q = xp2;
1657 }
1658 }
1659 wi[gas_comp_idx] = 1.0/(1.0/(trans_mult * this->well_index_[perf]) + (consistent_Q/2 * d / scaling));
1660
1661 return wi;
1662 }
1663
1664 template <typename TypeTag>
1665 void
1667 updateConnectionDFactor(const Simulator& simulator,
1668 SingleWellState<Scalar>& ws) const
1669 {
1670 if (! this->well_ecl_.getWDFAC().useDFactor()) {
1671 return;
1672 }
1673
1674 auto& d_factor = ws.perf_data.connection_d_factor;
1675
1676 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
1677 const int cell_idx = this->well_cells_[perf];
1678 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1679
1680 d_factor[perf] = this->computeConnectionDFactor(perf, intQuants, ws);
1681 }
1682 }
1683
1684 template <typename TypeTag>
1687 computeConnectionDFactor(const int perf,
1688 const IntensiveQuantities& intQuants,
1689 const SingleWellState<Scalar>& ws) const
1690 {
1691 auto rhoGS = [regIdx = this->pvtRegionIdx()]() {
1692 return FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regIdx);
1693 };
1694
1695 // Viscosity is evaluated at connection pressure.
1696 auto gas_visc = [connection_pressure = ws.perf_data.pressure[perf],
1697 temperature = ws.temperature,
1698 regIdx = this->pvtRegionIdx(), &intQuants]()
1699 {
1700 const auto rv = getValue(intQuants.fluidState().Rv());
1701
1702 const auto& gasPvt = FluidSystem::gasPvt();
1703
1704 // Note that rv here is from grid block with typically
1705 // p_block > connection_pressure
1706 // so we may very well have rv > rv_sat
1707 const Scalar rv_sat = gasPvt.saturatedOilVaporizationFactor
1708 (regIdx, temperature, connection_pressure);
1709
1710 if (! (rv < rv_sat)) {
1711 return gasPvt.saturatedViscosity(regIdx, temperature,
1712 connection_pressure);
1713 }
1714
1715 return gasPvt.viscosity(regIdx, temperature, connection_pressure,
1716 rv, getValue(intQuants.fluidState().Rvw()));
1717 };
1718
1719 const auto& connection = this->well_ecl_.getConnections()
1720 [ws.perf_data.ecl_index[perf]];
1721
1722 return this->well_ecl_.getWDFAC().getDFactor(rhoGS, gas_visc, connection);
1723 }
1724
1725
1726 template <typename TypeTag>
1727 void
1730 SingleWellState<Scalar>& ws) const
1731 {
1732 auto connCF = [&connIx = std::as_const(ws.perf_data.ecl_index),
1733 &conns = this->well_ecl_.getConnections()]
1734 (const int perf)
1735 {
1736 return conns[connIx[perf]].CF();
1737 };
1738
1739 auto& tmult = ws.perf_data.connection_compaction_tmult;
1740 auto& ctf = ws.perf_data.connection_transmissibility_factor;
1741
1742 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
1743 const int cell_idx = this->well_cells_[perf];
1744
1745 const auto& intQuants = simulator.model()
1746 .intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1747
1748 tmult[perf] = simulator.problem()
1749 .template wellTransMultiplier<double>(intQuants, cell_idx);
1750
1751 ctf[perf] = connCF(perf) * tmult[perf];
1752 }
1753 }
1754
1755
1756 template<typename TypeTag>
1759 {
1760 if constexpr (Indices::oilEnabled) {
1761 return fs.pressure(FluidSystem::oilPhaseIdx);
1762 } else if constexpr (Indices::gasEnabled) {
1763 return fs.pressure(FluidSystem::gasPhaseIdx);
1764 } else {
1765 return fs.pressure(FluidSystem::waterPhaseIdx);
1766 }
1767 }
1768
1769 template <typename TypeTag>
1770 template<class Value, class Callback>
1771 void
1773 getMobility(const Simulator& simulator,
1774 const int perf,
1775 std::vector<Value>& mob,
1776 Callback& extendEval,
1777 [[maybe_unused]] DeferredLogger& deferred_logger) const
1778 {
1779 auto relpermArray = []()
1780 {
1781 if constexpr (std::is_same_v<Value, Scalar>) {
1782 return std::array<Scalar,3>{};
1783 } else {
1784 return std::array<Eval,3>{};
1785 }
1786 };
1787 const int cell_idx = this->well_cells_[perf];
1788 assert (int(mob.size()) == this->num_components_);
1789 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
1790 const auto& materialLawManager = simulator.problem().materialLawManager();
1791
1792 // either use mobility of the perforation cell or calculate its own
1793 // based on passing the saturation table index
1794 const int satid = this->saturation_table_number_[perf] - 1;
1795 const int satid_elem = materialLawManager->satnumRegionIdx(cell_idx);
1796 if (satid == satid_elem) { // the same saturation number is used. i.e. just use the mobilty from the cell
1797 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
1798 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1799 continue;
1800 }
1801
1802 const unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
1803 mob[activeCompIdx] = extendEval(intQuants.mobility(phaseIdx));
1804 }
1805 if constexpr (has_solvent) {
1806 mob[Indices::contiSolventEqIdx] = extendEval(intQuants.solventMobility());
1807 }
1808 } else {
1809 const auto& paramsCell = materialLawManager->connectionMaterialLawParams(satid, cell_idx);
1810 auto relativePerms = relpermArray();
1811 MaterialLaw::relativePermeabilities(relativePerms, paramsCell, intQuants.fluidState());
1812
1813 // reset the satnumvalue back to original
1814 materialLawManager->connectionMaterialLawParams(satid_elem, cell_idx);
1815
1816 // compute the mobility
1817 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
1818 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1819 continue;
1820 }
1821
1822 const unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
1823 mob[activeCompIdx] = extendEval(relativePerms[phaseIdx] / intQuants.fluidState().viscosity(phaseIdx));
1824 }
1825
1826 // this may not work if viscosity and relperms has been modified?
1827 if constexpr (has_solvent) {
1828 OPM_DEFLOG_THROW(std::runtime_error, "individual mobility for wells does not work in combination with solvent", deferred_logger);
1829 }
1830 }
1831
1832 if (this->isInjector() && !this->inj_fc_multiplier_.empty()) {
1833 const auto perf_ecl_index = this->perforationData()[perf].ecl_index;
1834 const auto& connections = this->well_ecl_.getConnections();
1835 const auto& connection = connections[perf_ecl_index];
1836 if (connection.filterCakeActive()) {
1837 for (auto& val : mob) {
1838 val *= this->inj_fc_multiplier_[perf];
1839 }
1840 }
1841 }
1842 }
1843
1844
1845 template<typename TypeTag>
1846 bool
1849 WellState<Scalar>& well_state,
1850 DeferredLogger& deferred_logger) const
1851 {
1852 const auto& summary_state = simulator.vanguard().summaryState();
1853
1854 auto bhp_at_thp_limit = computeBhpAtThpLimitProdWithAlq(
1855 simulator, summary_state, this->getALQ(well_state), deferred_logger);
1856 if (bhp_at_thp_limit) {
1857 std::vector<Scalar> rates(this->number_of_phases_, 0.0);
1858 if (thp_update_iterations) {
1859 computeWellRatesWithBhpIterations(simulator, *bhp_at_thp_limit,
1860 rates, deferred_logger);
1861 } else {
1862 computeWellRatesWithBhp(simulator, *bhp_at_thp_limit,
1863 rates, deferred_logger);
1864 }
1865 auto& ws = well_state.well(this->name());
1866 ws.surface_rates = rates;
1867 ws.bhp = *bhp_at_thp_limit;
1868 ws.thp = this->getTHPConstraint(summary_state);
1869 return true;
1870 } else {
1871 return false;
1872 }
1873 }
1874
1875 template <typename TypeTag>
1876 void
1879 const std::function<Scalar(const Scalar)>& connPICalc,
1880 const std::vector<Scalar>& mobility,
1881 Scalar* connPI) const
1882 {
1883 const auto& pu = this->phaseUsage();
1884 const int np = this->number_of_phases_;
1885 for (int p = 0; p < np; ++p) {
1886 // Note: E100's notion of PI value phase mobility includes
1887 // the reciprocal FVF.
1888 const auto connMob =
1889 mobility[this->flowPhaseToModelCompIdx(p)]
1890 * fs.invB(this->flowPhaseToModelPhaseIdx(p)).value();
1891
1892 connPI[p] = connPICalc(connMob);
1893 }
1894
1895 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
1896 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
1897 {
1898 const auto io = pu.phase_pos[Oil];
1899 const auto ig = pu.phase_pos[Gas];
1900
1901 const auto vapoil = connPI[ig] * fs.Rv().value();
1902 const auto disgas = connPI[io] * fs.Rs().value();
1903
1904 connPI[io] += vapoil;
1905 connPI[ig] += disgas;
1906 }
1907 }
1908
1909
1910 template <typename TypeTag>
1911 void
1914 const Phase preferred_phase,
1915 const std::function<Scalar(const Scalar)>& connIICalc,
1916 const std::vector<Scalar>& mobility,
1917 Scalar* connII,
1918 DeferredLogger& deferred_logger) const
1919 {
1920 // Assumes single phase injection
1921 const auto& pu = this->phaseUsage();
1922
1923 auto phase_pos = 0;
1924 if (preferred_phase == Phase::GAS) {
1925 phase_pos = pu.phase_pos[Gas];
1926 }
1927 else if (preferred_phase == Phase::OIL) {
1928 phase_pos = pu.phase_pos[Oil];
1929 }
1930 else if (preferred_phase == Phase::WATER) {
1931 phase_pos = pu.phase_pos[Water];
1932 }
1933 else {
1934 OPM_DEFLOG_THROW(NotImplemented,
1935 fmt::format("Unsupported Injector Type ({}) "
1936 "for well {} during connection I.I. calculation",
1937 static_cast<int>(preferred_phase), this->name()),
1938 deferred_logger);
1939 }
1940
1941 const auto mt = std::accumulate(mobility.begin(), mobility.end(), 0.0);
1942 connII[phase_pos] = connIICalc(mt * fs.invB(this->flowPhaseToModelPhaseIdx(phase_pos)).value());
1943 }
1944
1945} // namespace Opm
#define OPM_DEFLOG_THROW(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:45
@ Liquid
Definition: BlackoilPhases.hpp:42
@ Aqua
Definition: BlackoilPhases.hpp:42
@ Vapour
Definition: BlackoilPhases.hpp:42
Definition: DeferredLogger.hpp:57
void info(const std::string &tag, const std::string &message)
void warning(const std::string &tag, const std::string &message)
void debug(const std::string &tag, const std::string &message)
Definition: GroupState.hpp:38
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
Definition: SingleWellState.hpp:42
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:110
Scalar temperature
Definition: SingleWellState.hpp:90
Class for computing BHP limits.
Definition: WellBhpThpCalculator.hpp:41
std::optional< Scalar > estimateStableBhp(const WellState< Scalar > &well_state, const Well &well, const std::vector< Scalar > &rates, const Scalar rho, const SummaryState &summaryState) const
Scalar calculateMinimumBhpFromThp(const WellState< Scalar > &well_state, const Well &well, const SummaryState &summaryState, const Scalar rho) const
bool isStableSolution(const WellState< Scalar > &well_state, const Well &well, const std::vector< Scalar > &rates, const SummaryState &summaryState) const
EvalWell calculateBhpFromThp(const WellState< Scalar > &well_state, const std::vector< EvalWell > &rates, const Well &well, const SummaryState &summaryState, const Scalar rho, DeferredLogger &deferred_logger) const
int number_of_perforations_
Definition: WellInterfaceGeneric.hpp:306
Well well_ecl_
Definition: WellInterfaceGeneric.hpp:273
const std::vector< FluidSystem::Scalar > & wellIndex() const
Definition: WellInterfaceGeneric.hpp:141
FluidSystem::Scalar wsolvent_
Definition: WellInterfaceGeneric.hpp:348
Definition: WellInterfaceIndices.hpp:34
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:195
Scalar woxygen() const
Definition: WellInterface_impl.hpp:165
IndividualOrGroup
Definition: WellInterface.hpp:244
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:78
WellInterface(const Well &well, const ParallelWellInfo< Scalar > &pw_info, const int time_step, const ModelParameters &param, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_components, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Constructor.
Definition: WellInterface_impl.hpp:55
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1773
void assembleWellEq(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:782
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:100
BlackOilFluidState< Eval, FluidSystem, has_temperature, has_energy, Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, Indices::numPhases > FluidState
Definition: WellInterface.hpp:132
Scalar wfoam() const
Definition: WellInterface_impl.hpp:127
std::vector< RateVector > connectionRates_
Definition: WellInterface.hpp:389
bool wellUnderZeroRateTarget(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1444
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:1878
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1848
void addCellRates(RateVector &rates, int cellIdx) const
Definition: WellInterface_impl.hpp:874
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:890
virtual void init(const PhaseUsage *phase_usage_arg, const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const int num_cells, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step)
Definition: WellInterface_impl.hpp:92
bool iterateWellEquations(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:481
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:82
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
void checkWellOperability(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:910
std::optional< Scalar > estimateOperableBhp(const Simulator &ebos_simulator, const double dt, WellState< Scalar > &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:611
void updateWellOperability(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1013
void updateWellStateRates(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1540
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupState< Scalar > &group_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1067
Eval getPerfCellPressure(const FluidState &fs) const
Definition: WellInterface_impl.hpp:1758
void updateConnectionTransmissibilityFactor(const Simulator &simulator, SingleWellState< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1729
bool stoppedOrZeroRateTarget(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1480
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:80
static constexpr bool has_solvent
Definition: WellInterface.hpp:108
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:816
bool gliftBeginTimeStepWellTestIterateWellEquations(const Simulator &ebos_simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:935
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:85
void computeConnLevelInjInd(const FluidState &fs, const Phase preferred_phase, const std::function< Scalar(const Scalar)> &connIICalc, const std::vector< Scalar > &mobility, Scalar *connII, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1913
void solveWellEquation(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:725
void updateConnectionDFactor(const Simulator &simulator, SingleWellState< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1667
std::vector< Scalar > initialWellRateFractions(const Simulator &ebosSimulator, const WellState< Scalar > &well_state) const
Definition: WellInterface_impl.hpp:1493
Scalar wsalt() const
Definition: WellInterface_impl.hpp:141
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1049
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:969
typename Base::Eval Eval
Definition: WellInterface.hpp:95
void wellTesting(const Simulator &simulator, const double simulation_time, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, WellTestState &welltest_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:376
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:111
Scalar computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellState< Scalar > &ws) const
Definition: WellInterface_impl.hpp:1687
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:634
bool solveWellWithZeroRate(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:670
bool solveWellForTesting(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:690
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:81
bool solveWellWithTHPConstraint(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:514
Scalar wurea() const
Definition: WellInterface_impl.hpp:183
bool wellUnderZeroGroupRateTarget(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger, std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1461
Scalar wmicrobes() const
Definition: WellInterface_impl.hpp:153
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:797
bool updateWellControlAndStatusLocalIteration(const Simulator &simulator, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const Scalar WQTotal, DeferredLogger &deferred_logger, const bool fixed_control=false, const bool fixed_status=false)
Definition: WellInterface_impl.hpp:271
static constexpr bool has_zFraction
Definition: WellInterface.hpp:109
Definition: WellState.hpp:62
const SingleWellState< Scalar > & well(std::size_t well_index) const
Definition: WellState.hpp:307
void setALQ(const std::string &name, Scalar value)
Definition: WellState.hpp:186
int numPhases() const
The number of phases present.
Definition: WellState.hpp:261
Scalar getALQ(const std::string &name) const
Definition: WellState.hpp:181
@ NONE
Definition: DeferredLogger.hpp:46
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:37
PhaseUsage phaseUsage(const Phases &phases)
Determine the active phases.
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:147
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: BlackoilPhases.hpp:46