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#ifndef OPM_WELLINTERFACE_IMPL_HEADER_INCLUDED
23#define OPM_WELLINTERFACE_IMPL_HEADER_INCLUDED
24
25// Improve IDE experience
26#ifndef OPM_WELLINTERFACE_HEADER_INCLUDED
27#include <config.h>
29#endif
30
31#include <opm/common/Exceptions.hpp>
32
33#include <opm/input/eclipse/Schedule/ScheduleTypes.hpp>
34#include <opm/input/eclipse/Schedule/Well/WDFAC.hpp>
35
37
42
43#include <dune/common/version.hh>
44
45#include <algorithm>
46#include <cassert>
47#include <cstddef>
48#include <numbers>
49#include <utility>
50
51#include <fmt/format.h>
52
53namespace Opm
54{
55
56
57 template<typename TypeTag>
59 WellInterface(const Well& well,
60 const ParallelWellInfo<Scalar>& pw_info,
61 const int time_step,
62 const ModelParameters& param,
63 const RateConverterType& rate_converter,
64 const int pvtRegionIdx,
65 const int num_conservation_quantities,
66 const int num_phases,
67 const int index_of_well,
68 const std::vector<PerforationData<Scalar>>& perf_data)
70 pw_info,
71 time_step,
72 param,
73 rate_converter,
74 pvtRegionIdx,
75 num_conservation_quantities,
76 num_phases,
77 index_of_well,
78 perf_data)
79 {
81
82 if constexpr (has_solvent || has_zFraction) {
83 if (well.isInjector()) {
84 auto injectorType = this->well_ecl_.injectorType();
85 if (injectorType == InjectorType::GAS) {
86 this->wsolvent_ = this->well_ecl_.getSolventFraction();
87 }
88 }
89 }
90 }
91
92
93 template<typename TypeTag>
94 void
96 init(const std::vector<Scalar>& /* depth_arg */,
97 const Scalar gravity_arg,
98 const std::vector<Scalar>& B_avg,
99 const bool changed_to_open_this_step)
100 {
101 this->gravity_ = gravity_arg;
102 B_avg_ = B_avg;
103 this->changed_to_open_this_step_ = changed_to_open_this_step;
104 }
105
106
107
108
109 template<typename TypeTag>
112 wpolymer() const
113 {
114 if constexpr (has_polymer) {
115 return this->wpolymer_();
116 }
117
118 return 0.0;
119 }
120
121
122
123
124
125 template<typename TypeTag>
128 wfoam() const
129 {
130 if constexpr (has_foam) {
131 return this->wfoam_();
132 }
133
134 return 0.0;
135 }
136
137
138
139 template<typename TypeTag>
142 wsalt() const
143 {
144 if constexpr (has_brine) {
145 return this->wsalt_();
146 }
147
148 return 0.0;
149 }
150
151 template<typename TypeTag>
154 wmicrobes() const
155 {
156 if constexpr (has_micp) {
157 return this->wmicrobes_();
158 }
159
160 return 0.0;
161 }
162
163 template<typename TypeTag>
166 woxygen() const
167 {
168 if constexpr (has_micp) {
169 return this->woxygen_();
170 }
171
172 return 0.0;
173 }
174
175 template<typename TypeTag>
178 wurea() const
179 {
180 if constexpr (has_micp) {
181 return this->wurea_();
182 }
183
184 return 0.0;
185 }
186
187 template<typename TypeTag>
188 bool
190 updateWellControl(const Simulator& simulator,
191 const IndividualOrGroup iog,
192 const GroupStateHelperType& groupStateHelper,
193 WellStateType& well_state) /* const */
194 {
195 auto& deferred_logger = groupStateHelper.deferredLogger();
196 OPM_TIMEFUNCTION();
197 if (stoppedOrZeroRateTarget(groupStateHelper)) {
198 return false;
199 }
200
201 const auto& summaryState = simulator.vanguard().summaryState();
202 const auto& schedule = simulator.vanguard().schedule();
203 const auto& well = this->well_ecl_;
204 auto& ws = well_state.well(this->index_of_well_);
205 std::string from;
206 bool is_grup = false;
207 if (well.isInjector()) {
208 from = WellInjectorCMode2String(ws.injection_cmode);
209 is_grup = ws.injection_cmode == Well::InjectorCMode::GRUP;
210 } else {
211 from = WellProducerCMode2String(ws.production_cmode);
212 is_grup = ws.production_cmode == Well::ProducerCMode::GRUP;
213 }
214
215 const int episodeIdx = simulator.episodeIndex();
216 const auto& iterCtx = simulator.problem().iterationContext();
217 const int nupcol = schedule[episodeIdx].nupcol();
218 const bool oscillating =
219 std::ranges::count(this->well_control_log_, from) >= this->param_.max_number_of_well_switches_;
220 if (oscillating && !is_grup) { // we would like to avoid ending up as GRUP
221 // only output first time
222 const bool output =
223 std::ranges::count(this->well_control_log_, from) == this->param_.max_number_of_well_switches_;
224 if (output) {
225 const auto msg = fmt::format(" The control mode for well {} is oscillating. \n"
226 "We don't allow for more than {} switches after NUPCOL iterations. (NUPCOL = {}) \n"
227 "The control is kept at {}.",
228 this->name(), this->param_.max_number_of_well_switches_, nupcol, from);
229 deferred_logger.info(msg);
230 // add one more to avoid outputting the same info again
231 this->well_control_log_.push_back(from);
232 }
233 return false;
234 }
235 bool changed = false;
236 if (iog == IndividualOrGroup::Individual) {
237 changed = this->checkIndividualConstraints(ws, summaryState, deferred_logger);
238 } else if (iog == IndividualOrGroup::Group) {
239 changed = this->checkGroupConstraints(
240 groupStateHelper, schedule, summaryState, true, well_state
241 );
242 } else {
243 assert(iog == IndividualOrGroup::Both);
244 changed = this->checkConstraints(groupStateHelper, schedule, summaryState, well_state);
245 }
246 Parallel::Communication cc = simulator.vanguard().grid().comm();
247 // checking whether control changed
248 if (changed) {
249 std::string to;
250 if (well.isInjector()) {
251 to = WellInjectorCMode2String(ws.injection_cmode);
252 } else {
253 to = WellProducerCMode2String(ws.production_cmode);
254 }
255 std::ostringstream ss;
256 ss << " Switching control mode for well " << this->name()
257 << " from " << from
258 << " to " << to;
259 if (iterCtx.inLocalSolve()) {
260 ss << " (NLDD domain solve)";
261 }
262 if (cc.size() > 1) {
263 ss << " on rank " << cc.rank();
264 }
265 deferred_logger.debug(ss.str());
266
267 // We always store the current control as it is used for output
268 // and only after iteration >= nupcol
269 // we log all switches to check if the well controls oscillates.
270 // Skip logging during NLDD local solves to avoid exhausting the
271 // global oscillation budget with provisional domain-level switches.
272 if (!iterCtx.inLocalSolve()) {
273 if (!iterCtx.withinNupcol(nupcol) || this->well_control_log_.empty()) {
274 this->well_control_log_.push_back(from);
275 }
276 }
277 updateWellStateWithTarget(simulator, groupStateHelper, well_state);
278 updatePrimaryVariables(groupStateHelper);
279 }
280
281 return changed;
282 }
283
284 template<typename TypeTag>
285 bool
288 const GroupStateHelperType& groupStateHelper,
289 const Well::InjectionControls& inj_controls,
290 const Well::ProductionControls& prod_controls,
291 const Scalar wqTotal,
292 WellStateType& well_state,
293 const bool fixed_control,
294 const bool fixed_status,
295 const bool solving_with_zero_rate)
296 {
297 OPM_TIMEFUNCTION();
298 auto& deferred_logger = groupStateHelper.deferredLogger();
299 const auto& summary_state = simulator.vanguard().summaryState();
300 const auto& schedule = simulator.vanguard().schedule();
301 auto& ws = well_state.well(this->index_of_well_);
302 std::string from;
303 if (this->isInjector()) {
304 from = WellInjectorCMode2String(ws.injection_cmode);
305 } else {
306 from = WellProducerCMode2String(ws.production_cmode);
307 }
308 const bool oscillating =
309 std::ranges::count(this->well_control_log_, from) >= this->param_.max_number_of_well_switches_;
310
311 if (oscillating || this->wellUnderZeroRateTarget(groupStateHelper) || !(well_state.well(this->index_of_well_).status == WellStatus::OPEN)) {
312 return false;
313 }
314
315 const Scalar sgn = this->isInjector() ? 1.0 : -1.0;
316 if (!this->wellIsStopped()){
317 if (wqTotal*sgn <= 0.0 && !fixed_status){
318 this->stopWell();
319 return true;
320 } else {
321 bool changed = false;
322 if (!fixed_control) {
323 // When solving_with_zero_rate=true, fixed_control=true, so this block should never
324 // be entered.
325 if (solving_with_zero_rate) {
326 OPM_DEFLOG_THROW(std::runtime_error, fmt::format(
327 "Well {}: solving_with_zero_rate should not be true when fixed_control is false",
328 this->name()), deferred_logger);
329 }
330
331 // Changing to group controls here may lead to inconsistencies in the group handling which in turn
332 // may result in excessive back and forth switching. However, we currently allow this by default.
333 // The switch check_group_constraints_inner_well_iterations_ is a temporary solution.
334 const bool hasGroupControl = this->isInjector() ? inj_controls.hasControl(Well::InjectorCMode::GRUP) :
335 prod_controls.hasControl(Well::ProducerCMode::GRUP);
336 bool isGroupControl = ws.production_cmode == Well::ProducerCMode::GRUP || ws.injection_cmode == Well::InjectorCMode::GRUP;
337 if (! (isGroupControl && !this->param_.check_group_constraints_inner_well_iterations_)) {
338 changed = this->checkIndividualConstraints(ws, summary_state, deferred_logger, inj_controls, prod_controls);
339 }
340 if (hasGroupControl && this->param_.check_group_constraints_inner_well_iterations_) {
341 changed = changed || this->checkGroupConstraints(
342 groupStateHelper, schedule, summary_state, false, well_state
343 );
344 }
345
346 if (changed) {
347 const bool thp_controlled = this->isInjector() ? ws.injection_cmode == Well::InjectorCMode::THP :
348 ws.production_cmode == Well::ProducerCMode::THP;
349 if (thp_controlled){
350 ws.thp = this->getTHPConstraint(summary_state);
351 } else {
352 // don't call for thp since this might trigger additional local solve
353 updateWellStateWithTarget(simulator, groupStateHelper, well_state);
354 }
355 updatePrimaryVariables(groupStateHelper);
356 }
357 }
358 return changed;
359 }
360 } else if (!fixed_status){
361 // well is stopped, check if current bhp allows reopening
362 const Scalar bhp = well_state.well(this->index_of_well_).bhp;
363 Scalar prod_limit = prod_controls.bhp_limit;
364 Scalar inj_limit = inj_controls.bhp_limit;
365 const bool has_thp = this->wellHasTHPConstraints(summary_state);
366 if (has_thp){
367 std::vector<Scalar> rates(this->num_conservation_quantities_);
368 if (this->isInjector()){
369 const Scalar bhp_thp = WellBhpThpCalculator(*this).
370 calculateBhpFromThp(well_state, rates,
371 this->well_ecl_,
372 summary_state,
373 this->getRefDensity(),
374 deferred_logger);
375 inj_limit = std::min(bhp_thp, static_cast<Scalar>(inj_controls.bhp_limit));
376 } else {
377 // if the well can operate, it must at least be able to produce
378 // at the lowest bhp of the bhp-curve (explicit fractions)
379 const Scalar bhp_min = WellBhpThpCalculator(*this).
380 calculateMinimumBhpFromThp(well_state,
381 this->well_ecl_,
382 summary_state,
383 this->getRefDensity());
384 prod_limit = std::max(bhp_min, static_cast<Scalar>(prod_controls.bhp_limit));
385 }
386 }
387 const Scalar bhp_diff = (this->isInjector())? inj_limit - bhp: bhp - prod_limit;
388 if (bhp_diff > 0){
389 this->openWell();
390 well_state.well(this->index_of_well_).bhp = (this->isInjector())? inj_limit : prod_limit;
391 if (has_thp) {
392 well_state.well(this->index_of_well_).thp = this->getTHPConstraint(summary_state);
393 }
394 return true;
395 } else {
396 return false;
397 }
398 } else {
399 return false;
400 }
401 }
402
403 template<typename TypeTag>
404 void
406 wellTesting(const Simulator& simulator,
407 const double simulation_time,
408 const GroupStateHelperType& groupStateHelper,
409 WellStateType& well_state,
410 WellTestState& well_test_state,
411 GLiftEclWells& ecl_well_map,
412 std::map<std::string, double>& open_times)
413 {
414 OPM_TIMEFUNCTION();
415 auto& deferred_logger = groupStateHelper.deferredLogger();
416 const auto& group_state = groupStateHelper.groupState();
417 deferred_logger.info(" well " + this->name() + " is being tested");
418
419 GroupStateHelperType groupStateHelper_copy = groupStateHelper;
420 WellStateType well_state_copy = well_state;
421 // Ensure that groupStateHelper uses well_state_copy as WellState for the well testing
422 // and the guard ensures that the original well state is restored at scope exit, i.e. at
423 // the end of this function.
424 auto guard = groupStateHelper_copy.pushWellState(well_state_copy);
425 auto& ws = well_state_copy.well(this->indexOfWell());
426
427 const auto& summary_state = simulator.vanguard().summaryState();
428 const bool has_thp_limit = this->wellHasTHPConstraints(summary_state);
429 if (this->isProducer()) {
430 ws.production_cmode = has_thp_limit ? Well::ProducerCMode::THP : Well::ProducerCMode::BHP;
431 } else {
432 ws.injection_cmode = has_thp_limit ? Well::InjectorCMode::THP : Well::InjectorCMode::BHP;
433 }
434 // We test the well as an open well during the well testing
435 ws.open();
436
437 scaleSegmentRatesAndPressure(well_state_copy);
438 calculateExplicitQuantities(simulator, groupStateHelper_copy);
439 updatePrimaryVariables(groupStateHelper_copy);
440
441 if (this->isProducer()) {
442 const auto& schedule = simulator.vanguard().schedule();
443 const auto report_step = simulator.episodeIndex();
444 const auto& glo = schedule.glo(report_step);
445 if (glo.active()) {
446 gliftBeginTimeStepWellTestUpdateALQ(simulator,
447 well_state_copy,
448 group_state,
449 ecl_well_map,
450 deferred_logger);
451 }
452 }
453
454 WellTestState welltest_state_temp;
455
456 bool testWell = true;
457 // if a well is closed because all completions are closed, we need to check each completion
458 // individually. We first open all completions, then we close one by one by calling updateWellTestState
459 // untill the number of closed completions do not increase anymore.
460 while (testWell) {
461 const std::size_t original_number_closed_completions = welltest_state_temp.num_closed_completions();
462 bool converged = solveWellForTesting(simulator, groupStateHelper_copy, well_state_copy);
463 if (!converged) {
464 const auto msg = fmt::format("WTEST: Well {} is not solvable (physical)", this->name());
465 deferred_logger.debug(msg);
466 return;
467 }
468
469
470 updateWellOperability(simulator, well_state_copy, groupStateHelper_copy);
471 if ( !this->isOperableAndSolvable() ) {
472 const auto msg = fmt::format("WTEST: Well {} is not operable (physical)", this->name());
473 deferred_logger.debug(msg);
474 return;
475 }
476 std::vector<Scalar> potentials;
477 try {
478 computeWellPotentials(simulator, well_state_copy, groupStateHelper_copy, potentials);
479 } catch (const std::exception& e) {
480 const std::string msg = fmt::format("well {}: computeWellPotentials() "
481 "failed during testing for re-opening: ",
482 this->name(), e.what());
483 deferred_logger.info(msg);
484 return;
485 }
486 const int np = well_state_copy.numPhases();
487 for (int p = 0; p < np; ++p) {
488 ws.well_potentials[p] = std::max(Scalar{0.0}, potentials[p]);
489 }
490 const bool under_zero_target = this->wellUnderZeroGroupRateTarget(groupStateHelper_copy);
491 this->updateWellTestState(well_state_copy.well(this->indexOfWell()),
492 simulation_time,
493 /*writeMessageToOPMLog=*/ false,
494 under_zero_target,
495 welltest_state_temp,
496 deferred_logger);
497 this->closeCompletions(welltest_state_temp);
498
499 // Stop testing if the well is closed or shut due to all completions shut
500 // Also check if number of completions has increased. If the number of closed completions do not increased
501 // we stop the testing.
502 // TODO: it can be tricky here, if the well is shut/closed due to other reasons
503 if ( welltest_state_temp.num_closed_wells() > 0 ||
504 (original_number_closed_completions == welltest_state_temp.num_closed_completions()) ) {
505 testWell = false; // this terminates the while loop
506 }
507 }
508
509 // update wellTestState if the well test succeeds
510 if (!welltest_state_temp.well_is_closed(this->name())) {
511 well_test_state.open_well(this->name());
512
513 std::string msg = std::string("well ") + this->name() + std::string(" is re-opened");
514 deferred_logger.info(msg);
515
516 // also reopen completions
517 for (const auto& completion : this->well_ecl_.getCompletions()) {
518 if (!welltest_state_temp.completion_is_closed(this->name(), completion.first))
519 well_test_state.open_completion(this->name(), completion.first);
520 }
521 well_state = well_state_copy;
522 open_times.try_emplace(this->name(), well_test_state.lastTestTime(this->name()));
523 }
524 }
525
526
527
528
529 template<typename TypeTag>
530 bool
532 iterateWellEquations(const Simulator& simulator,
533 const double dt,
534 const GroupStateHelperType& groupStateHelper,
535 WellStateType& well_state)
536 {
537 OPM_TIMEFUNCTION();
538 auto& deferred_logger = groupStateHelper.deferredLogger();
539
540 const auto& summary_state = simulator.vanguard().summaryState();
541 const auto inj_controls = this->well_ecl_.isInjector() ? this->well_ecl_.injectionControls(summary_state) : Well::InjectionControls(0);
542 const auto prod_controls = this->well_ecl_.isProducer() ? this->well_ecl_.productionControls(summary_state) : Well::ProductionControls(0);
543 const auto& ws = well_state.well(this->indexOfWell());
544 const auto pmode_orig = ws.production_cmode;
545 const auto imode_orig = ws.injection_cmode;
546 bool converged = false;
547 try {
548 // TODO: the following two functions will be refactored to be one to reduce the code duplication
549 if (!this->param_.local_well_solver_control_switching_){
550 converged = this->iterateWellEqWithControl(simulator, dt, inj_controls, prod_controls, groupStateHelper, well_state);
551 } else {
552 if (this->param_.use_implicit_ipr_ && this->well_ecl_.isProducer() && (well_state.well(this->index_of_well_).status == WellStatus::OPEN)) {
553 converged = solveWellWithOperabilityCheck(
554 simulator, dt, inj_controls, prod_controls, groupStateHelper, well_state
555 );
556 } else {
557 converged = this->iterateWellEqWithSwitching(
558 simulator, dt, inj_controls, prod_controls, groupStateHelper, well_state,
559 /*fixed_control=*/false, /*fixed_status=*/false, /*solving_with_zero_rate=*/false
560 );
561 }
562 }
563
564 } catch (NumericalProblem& e ) {
565 const std::string msg = "Inner well iterations failed for well " + this->name() + " Treat the well as unconverged. ";
566 deferred_logger.warning("INNER_ITERATION_FAILED", msg);
567 converged = false;
568 }
569 if (converged) {
570 // Add debug info for switched controls
571 if (ws.production_cmode != pmode_orig || ws.injection_cmode != imode_orig) {
572 std::string from,to;
573 if (this->isInjector()) {
574 from = WellInjectorCMode2String(imode_orig);
575 to = WellInjectorCMode2String(ws.injection_cmode);
576 } else {
577 from = WellProducerCMode2String(pmode_orig);
578 to = WellProducerCMode2String(ws.production_cmode);
579 }
580 const auto msg = fmt::format(" Well {} switched from {} to {} during local solve", this->name(), from, to);
581 deferred_logger.debug(msg);
582 const int episodeIdx = simulator.episodeIndex();
583 const auto& iterCtx = simulator.problem().iterationContext();
584 const auto& schedule = simulator.vanguard().schedule();
585 const int nupcol = schedule[episodeIdx].nupcol();
586 // We always store the current control as it is used for output
587 // and only after iteration >= nupcol
588 // we log all switches to check if the well controls oscillates
589 if (!iterCtx.withinNupcol(nupcol) || this->well_control_log_.empty()) {
590 this->well_control_log_.push_back(from);
591 }
592 }
593 }
594
595 return converged;
596 }
597
598 template<typename TypeTag>
599 bool
602 const double dt,
603 const Well::InjectionControls& inj_controls,
604 const Well::ProductionControls& prod_controls,
605 const GroupStateHelperType& groupStateHelper,
606 WellStateType& well_state)
607 {
608 OPM_TIMEFUNCTION();
609 auto& deferred_logger = groupStateHelper.deferredLogger();
610
611 const auto& summary_state = simulator.vanguard().summaryState();
612 bool converged = true;
613 auto& ws = well_state.well(this->index_of_well_);
614 // if well is stopped, check if we can reopen with explicit fraction
615 if (this->wellIsStopped()) {
616 this->openWell();
617 const bool use_vfpexplicit = this->operability_status_.use_vfpexplicit;
618 this->operability_status_.use_vfpexplicit = true;
619 auto bhp_target = estimateOperableBhp(simulator, dt, groupStateHelper, summary_state, well_state);
620 if (!bhp_target.has_value()) {
621 // no intersection with ipr
622 const auto msg = fmt::format("estimateOperableBhp: Did not find operable BHP for well {}", this->name());
623 deferred_logger.debug(msg);
624 // well can't operate using explicit fractions stop the well
625 // solve with zero rates
626 converged = solveWellWithZeroRate(simulator, dt, groupStateHelper, well_state);
627 this->stopWell();
628 this->operability_status_.can_obtain_bhp_with_thp_limit = false;
629 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
630 return converged;
631 } else {
632 // solve well with the estimated target bhp (or limit)
633 ws.thp = this->getTHPConstraint(summary_state);
634 const Scalar bhp = std::max(bhp_target.value(),
635 static_cast<Scalar>(prod_controls.bhp_limit));
636 solveWellWithBhp(simulator, dt, bhp, groupStateHelper, well_state);
637 this->operability_status_.use_vfpexplicit = use_vfpexplicit;
638 }
639 }
640 // solve well-equation
641 converged = this->iterateWellEqWithSwitching(
642 simulator, dt, inj_controls, prod_controls, groupStateHelper, well_state,
643 /*fixed_control=*/false, /*fixed_status=*/false, /*solving_with_zero_rate=*/false
644 );
645
646
647 const bool isThp = ws.production_cmode == Well::ProducerCMode::THP;
648 // check stability of solution under thp-control
649 if (converged && !stoppedOrZeroRateTarget(groupStateHelper) && isThp) {
650 auto rates = well_state.well(this->index_of_well_).surface_rates;
651 this->adaptRatesForVFP(rates);
652 this->updateIPRImplicit(simulator, groupStateHelper, well_state);
653 bool is_stable = WellBhpThpCalculator(*this).isStableSolution(well_state, this->well_ecl_, rates, summary_state);
654 if (!is_stable) {
655 // solution converged to an unstable point!
656 this->operability_status_.use_vfpexplicit = true;
657 auto bhp_stable = WellBhpThpCalculator(*this).estimateStableBhp(well_state, this->well_ecl_, rates, this->getRefDensity(), summary_state);
658 // if we find an intersection with a sufficiently lower bhp, re-solve equations
659 const Scalar reltol = 1e-3;
660 const Scalar cur_bhp = ws.bhp;
661 if (bhp_stable.has_value() && cur_bhp - bhp_stable.value() > cur_bhp*reltol){
662 const auto msg = fmt::format("Well {} converged to an unstable solution, re-solving", this->name());
663 deferred_logger.debug(msg);
664 solveWellWithBhp(
665 simulator, dt, bhp_stable.value(), groupStateHelper, well_state
666 );
667 // re-solve with hopefully good initial guess
668 ws.thp = this->getTHPConstraint(summary_state);
669 converged = this->iterateWellEqWithSwitching(
670 simulator, dt, inj_controls, prod_controls, groupStateHelper, well_state,
671 /*fixed_control=*/false, /*fixed_status=*/false, /*solving_with_zero_rate=*/false
672 );
673 }
674 }
675 }
676
677 if (!converged) {
678 // Well did not converge, switch to explicit fractions
679 this->operability_status_.use_vfpexplicit = true;
680 this->openWell();
681 auto bhp_target = estimateOperableBhp(
682 simulator, dt, groupStateHelper, summary_state, well_state
683 );
684 if (!bhp_target.has_value()) {
685 // solve with zero rate
686 // well can't operate using explicit fractions stop the well
687 converged = solveWellWithZeroRate(simulator, dt, groupStateHelper, well_state);
688 this->stopWell();
689 this->operability_status_.can_obtain_bhp_with_thp_limit = false;
690 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
691 return converged;
692 } else {
693 // solve well with the estimated target bhp (or limit)
694 const Scalar bhp = std::max(bhp_target.value(),
695 static_cast<Scalar>(prod_controls.bhp_limit));
696 solveWellWithBhp(
697 simulator, dt, bhp, groupStateHelper, well_state
698 );
699 ws.thp = this->getTHPConstraint(summary_state);
700 const auto msg = fmt::format("Well {} did not converge, re-solving with explicit fractions for VFP caculations.", this->name());
701 deferred_logger.debug(msg);
702 converged = this->iterateWellEqWithSwitching(simulator, dt,
703 inj_controls,
704 prod_controls,
705 groupStateHelper,
706 well_state,
707 /*fixed_control=*/false,
708 /*fixed_status=*/false,
709 /*solving_with_zero_rate=*/false);
710 }
711 }
712 // update operability
713 this->operability_status_.can_obtain_bhp_with_thp_limit = !this->wellIsStopped();
714 this->operability_status_.obey_thp_limit_under_bhp_limit = !this->wellIsStopped();
715 return converged;
716 }
717
718 template<typename TypeTag>
719 std::optional<typename WellInterface<TypeTag>::Scalar>
721 estimateOperableBhp(const Simulator& simulator,
722 const double dt,
723 const GroupStateHelperType& groupStateHelper,
724 const SummaryState& summary_state,
725 WellStateType& well_state)
726 {
727 if (!this->wellHasTHPConstraints(summary_state)) {
728 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summary_state);
729 const bool converged = solveWellWithBhp(
730 simulator, dt, bhp_limit, groupStateHelper, well_state
731 );
732 if (!converged || this->wellIsStopped()) {
733 return std::nullopt;
734 }
735
736 return bhp_limit;
737 }
738 OPM_TIMEFUNCTION();
739 // Given an unconverged well or closed well, estimate an operable bhp (if any)
740 // Get minimal bhp from vfp-curve
741 Scalar bhp_min = WellBhpThpCalculator(*this).calculateMinimumBhpFromThp(well_state, this->well_ecl_, summary_state, this->getRefDensity());
742 // Solve
743 const bool converged = solveWellWithBhp(
744 simulator, dt, bhp_min, groupStateHelper, well_state
745 );
746 if (!converged || this->wellIsStopped()) {
747 return std::nullopt;
748 }
749 this->updateIPRImplicit(simulator, groupStateHelper, well_state);
750 auto rates = well_state.well(this->index_of_well_).surface_rates;
751 this->adaptRatesForVFP(rates);
752 return WellBhpThpCalculator(*this).estimateStableBhp(well_state, this->well_ecl_, rates, this->getRefDensity(), summary_state);
753 }
754
755 template<typename TypeTag>
756 bool
758 solveWellWithBhp(const Simulator& simulator,
759 const double dt,
760 const Scalar bhp,
761 const GroupStateHelperType& groupStateHelper,
762 WellStateType& well_state)
763 {
764 OPM_TIMEFUNCTION();
765
766 // Solve a well using single bhp-constraint (but close if not operable under this)
767 auto group_state = GroupState<Scalar>(); // empty group
768 GroupStateHelperType groupStateHelper_copy = groupStateHelper;
769 // Ensure that groupStateHelper_copy uses the empty group state as GroupState for iterateWellEqWithSwitching()
770 // and the guard ensures that the original group state is restored at scope exit, i.e. at
771 // the end of this function.
772 auto group_guard = groupStateHelper_copy.pushGroupState(group_state);
773
774 auto inj_controls = Well::InjectionControls(0);
775 auto prod_controls = Well::ProductionControls(0);
776 auto& ws = well_state.well(this->index_of_well_);
777 auto cmode_inj = ws.injection_cmode;
778 auto cmode_prod = ws.production_cmode;
779 if (this->isInjector()) {
780 inj_controls.addControl(Well::InjectorCMode::BHP);
781 inj_controls.bhp_limit = bhp;
782 inj_controls.cmode = Well::InjectorCMode::BHP;
783 ws.injection_cmode = Well::InjectorCMode::BHP;
784 } else {
785 prod_controls.addControl(Well::ProducerCMode::BHP);
786 prod_controls.bhp_limit = bhp;
787 prod_controls.cmode = Well::ProducerCMode::BHP;
788 ws.production_cmode = Well::ProducerCMode::BHP;
789 }
790 // update well-state
791 ws.bhp = bhp;
792 // solve
793 const bool converged = this->iterateWellEqWithSwitching(
794 simulator, dt, inj_controls, prod_controls, groupStateHelper_copy,
795 well_state,
796 /*fixed_control=*/true,
797 /*fixed_status=*/false,
798 /*solving_with_zero_rate=*/false
799 );
800 ws.injection_cmode = cmode_inj;
801 ws.production_cmode = cmode_prod;
802 return converged;
803 }
804
805 template<typename TypeTag>
806 bool
808 solveWellWithZeroRate(const Simulator& simulator,
809 const double dt,
810 const GroupStateHelperType& groupStateHelper,
811 WellStateType& well_state)
812 {
813 OPM_TIMEFUNCTION();
814
815 // Solve a well as stopped with isolation (empty group state for assembly)
816 const auto well_status_orig = this->wellStatus_;
817 this->stopWell();
818
819 auto inj_controls = Well::InjectionControls(0);
820 auto prod_controls = Well::ProductionControls(0);
821
822 // Solve with well isolation - the flag "solving_with_zero_rate=true" will be passed down to
823 // assembleWellEqWithoutIterationImpl() to create an empty group state when assembling the
824 // well equations.
825 const bool converged = this->iterateWellEqWithSwitching(
826 simulator, dt, inj_controls, prod_controls,
827 groupStateHelper,
828 well_state,
829 /*fixed_control*/true,
830 /*fixed_status*/true,
831 /*solving_with_zero_rate*/true
832 );
833 this->wellStatus_ = well_status_orig;
834 return converged;
835 }
836
837 template<typename TypeTag>
838 bool
840 solveWellForTesting(const Simulator& simulator,
841 const GroupStateHelperType& groupStateHelper,
842 WellStateType& well_state)
843 {
844 OPM_TIMEFUNCTION();
845 auto& deferred_logger = groupStateHelper.deferredLogger();
846
847 const double dt = simulator.timeStepSize();
848
849 const auto& summary_state = simulator.vanguard().summaryState();
850 auto inj_controls = this->well_ecl_.isInjector() ? this->well_ecl_.injectionControls(summary_state) : Well::InjectionControls(0);
851 auto prod_controls = this->well_ecl_.isProducer() ? this->well_ecl_.productionControls(summary_state) : Well::ProductionControls(0);
852 this->onlyKeepBHPandTHPcontrols(summary_state, well_state, inj_controls, prod_controls);
853
854 bool converged = false;
855 try {
856 // TODO: the following two functions will be refactored to be one to reduce the code duplication
857 if (!this->param_.local_well_solver_control_switching_){
858 converged = this->iterateWellEqWithControl(
859 simulator, dt, inj_controls, prod_controls, groupStateHelper, well_state
860 );
861 } else {
862 if (this->param_.use_implicit_ipr_ && this->well_ecl_.isProducer() && (well_state.well(this->index_of_well_).status == WellStatus::OPEN)) {
863 converged = this->solveWellWithOperabilityCheck(
864 simulator, dt, inj_controls, prod_controls, groupStateHelper, well_state
865 );
866 } else {
867 converged = this->iterateWellEqWithSwitching(
868 simulator, dt, inj_controls, prod_controls, groupStateHelper, well_state,
869 /*fixed_control=*/false,
870 /*fixed_status=*/false,
871 /*solving_with_zero_rate=*/false
872 );
873 }
874 }
875
876 } catch (NumericalProblem& e ) {
877 const std::string msg = "Inner well iterations failed for well " + this->name() + " Treat the well as unconverged. ";
878 deferred_logger.warning("INNER_ITERATION_FAILED", msg);
879 converged = false;
880 }
881
882 if (converged) {
883 deferred_logger.debug("WellTest: Well equation for well " + this->name() + " converged");
884 return true;
885 }
886 const int max_iter = this->param_.max_welleq_iter_;
887 deferred_logger.debug("WellTest: Well equation for well " + this->name() + " failed converging in "
888 + std::to_string(max_iter) + " iterations");
889 return false;
890 }
891
892
893 template<typename TypeTag>
894 void
896 solveWellEquation(const Simulator& simulator,
897 const GroupStateHelperType& groupStateHelper,
898 WellStateType& well_state)
899 {
900 OPM_TIMEFUNCTION();
901 auto& deferred_logger = groupStateHelper.deferredLogger();
902 if (!this->isOperableAndSolvable() && !this->wellIsStopped())
903 return;
904
905 // keep a copy of the original well state
906 const WellStateType well_state0 = well_state;
907 const double dt = simulator.timeStepSize();
908 bool converged = iterateWellEquations(simulator, dt, groupStateHelper, well_state);
909
910 // Newly opened wells with THP control sometimes struggles to
911 // converge due to bad initial guess. Or due to the simple fact
912 // that the well needs to change to another control.
913 // We therefore try to solve the well with BHP control to get
914 // an better initial guess.
915 // If the well is supposed to operate under THP control
916 // "updateWellControl" will switch it back to THP later.
917 if (!converged) {
918 auto& ws = well_state.well(this->indexOfWell());
919 bool thp_control = false;
920 if (this->well_ecl_.isInjector()) {
921 thp_control = ws.injection_cmode == Well::InjectorCMode::THP;
922 if (thp_control) {
923 ws.injection_cmode = Well::InjectorCMode::BHP;
924 if (this->well_control_log_.empty()) { // only log the first control
925 this->well_control_log_.push_back(WellInjectorCMode2String(Well::InjectorCMode::THP));
926 }
927 }
928 } else {
929 thp_control = ws.production_cmode == Well::ProducerCMode::THP;
930 if (thp_control) {
931 ws.production_cmode = Well::ProducerCMode::BHP;
932 if (this->well_control_log_.empty()) { // only log the first control
933 this->well_control_log_.push_back(WellProducerCMode2String(Well::ProducerCMode::THP));
934 }
935 }
936 }
937 if (thp_control) {
938 const std::string msg = std::string("The newly opened well ") + this->name()
939 + std::string(" with THP control did not converge during inner iterations, we try again with bhp control");
940 deferred_logger.debug(msg);
941 converged = this->iterateWellEquations(simulator, dt, groupStateHelper, well_state);
942 }
943 }
944
945 if (!converged) {
946 const int max_iter = this->param_.max_welleq_iter_;
947 deferred_logger.debug("Compute initial well solution for well " + this->name() + ". Failed to converge in "
948 + std::to_string(max_iter) + " iterations");
949 well_state = well_state0;
950 }
951 }
952
953
954
955 template <typename TypeTag>
956 void
958 assembleWellEq(const Simulator& simulator,
959 const double dt,
960 const GroupStateHelperType& groupStateHelper,
961 WellStateType& well_state)
962 {
963 OPM_TIMEFUNCTION();
964 prepareWellBeforeAssembling(simulator, dt, groupStateHelper, well_state);
965 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, well_state,
966 /*solving_with_zero_rate=*/false);
967 }
968
969
970
971 template <typename TypeTag>
972 void
975 const GroupStateHelperType& groupStateHelper,
976 const double dt,
977 WellStateType& well_state,
978 const bool solving_with_zero_rate)
979 {
980 OPM_TIMEFUNCTION();
981 const auto& summary_state = simulator.vanguard().summaryState();
982 const auto inj_controls = this->well_ecl_.isInjector() ? this->well_ecl_.injectionControls(summary_state) : Well::InjectionControls(0);
983 const auto prod_controls = this->well_ecl_.isProducer() ? this->well_ecl_.productionControls(summary_state) : Well::ProductionControls(0);
984 // 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
985 // TODO: maybe we can use std::optional or pointers to simplify here
986 assembleWellEqWithoutIteration(simulator, groupStateHelper, dt, inj_controls, prod_controls, well_state, solving_with_zero_rate);
987 }
988
989
990
991 template<typename TypeTag>
992 void
995 const double dt,
996 const GroupStateHelperType& groupStateHelper,
997 WellStateType& well_state)
998 {
999 OPM_TIMEFUNCTION();
1000 auto& deferred_logger = groupStateHelper.deferredLogger();
1001 const bool old_well_operable = this->operability_status_.isOperableAndSolvable();
1002
1003 if (this->param_.check_well_operability_iter_)
1004 checkWellOperability(simulator, well_state, groupStateHelper);
1005
1006 // only use inner well iterations for the first newton iterations.
1007 const auto& iterCtx = simulator.problem().iterationContext();
1008 if (iterCtx.shouldRunInnerWellIterations(this->param_.max_niter_inner_well_iter_)) {
1009 const auto& ws = well_state.well(this->indexOfWell());
1010 const bool nonzero_rate_original =
1011 std::any_of(ws.surface_rates.begin(),
1012 ws.surface_rates.begin() + well_state.numPhases(),
1013 [](Scalar rate) { return rate != Scalar(0.0); });
1014
1015 this->operability_status_.solvable = true;
1016 if (number_of_well_reopenings_ >= this->param_.max_well_status_switch_) {
1017 // only output the first time
1018 if (number_of_well_reopenings_ == this->param_.max_well_status_switch_) {
1019 const std::string msg = fmt::format("well {} is oscillating between open and stop. \n"
1020 "We don't allow for more than {} re-openings "
1021 "and the well is therefore kept stopped.",
1022 this->name(), number_of_well_reopenings_);
1023 deferred_logger.debug(msg);
1024 changed_to_stopped_this_step_ = old_well_operable;
1025 } else {
1026 changed_to_stopped_this_step_ = false;
1027 }
1028 this->stopWell();
1029 bool converged_zero_rate = this->solveWellWithZeroRate(
1030 simulator, dt, groupStateHelper, well_state
1031 );
1032 if (this->param_.shut_unsolvable_wells_ && !converged_zero_rate ) {
1033 this->operability_status_.solvable = false;
1034 } else {
1035 this->operability_status_.can_obtain_bhp_with_thp_limit = false;
1036 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
1037 }
1038 // we increse the number of reopenings to avoid output in the next iteration
1039 number_of_well_reopenings_++;
1040 return;
1041 }
1042 bool converged = this->iterateWellEquations(
1043 simulator, dt, groupStateHelper, well_state
1044 );
1045
1046 if (converged) {
1047 const bool zero_target = this->wellUnderZeroRateTarget(groupStateHelper);
1048 if (this->wellIsStopped() && !zero_target && nonzero_rate_original) {
1049 // Well had non-zero rate, but was stopped during local well-solve. We re-open the well
1050 // for the next global iteration, but if the zero rate persists, it will be stopped.
1051 // This logic is introduced to prevent/ameliorate stopped/revived oscillations
1052 this->operability_status_.resetOperability();
1053 this->openWell();
1054 deferred_logger.debug(" " + this->name() + " is re-opened after being stopped during local solve");
1055 number_of_well_reopenings_++;
1056 }
1057 } else {
1058 // unsolvable wells are treated as not operable and will not be solved for in this iteration.
1059 if (this->param_.shut_unsolvable_wells_) {
1060 this->operability_status_.solvable = false;
1061 }
1062 }
1063 }
1064 if (this->operability_status_.has_negative_potentials) {
1065 auto well_state_copy = well_state;
1066 std::vector<Scalar> potentials;
1067 try {
1068 computeWellPotentials(simulator, well_state_copy, groupStateHelper, potentials);
1069 } catch (const std::exception& e) {
1070 const std::string msg = fmt::format("well {}: computeWellPotentials() failed "
1071 "during attempt to recompute potentials for well: ",
1072 this->name(), e.what());
1073 deferred_logger.info(msg);
1074 this->operability_status_.has_negative_potentials = true;
1075 }
1076 auto& ws = well_state.well(this->indexOfWell());
1077 const int np = well_state.numPhases();
1078 for (int p = 0; p < np; ++p) {
1079 ws.well_potentials[p] = std::max(Scalar{0.0}, potentials[p]);
1080 }
1081 }
1082 this->changed_to_open_this_step_ = false;
1083 changed_to_stopped_this_step_ = false;
1084
1085 const bool well_operable = this->operability_status_.isOperableAndSolvable();
1086 if (!well_operable) {
1087 this->stopWell();
1088 try {
1089 this->solveWellWithZeroRate(
1090 simulator, dt, groupStateHelper, well_state
1091 );
1092 } catch (const std::exception& e) {
1093 const std::string msg = fmt::format("well {}: solveWellWithZeroRate() failed "
1094 "during attempt to solve with zero rate for well: ",
1095 this->name(), e.what());
1096 deferred_logger.info(msg);
1097 // we set the rate to zero to make sure the well dont contribute to the group rate
1098 auto& ws = well_state.well(this->indexOfWell());
1099 const int np = well_state.numPhases();
1100 for (int p = 0; p < np; ++p) {
1101 ws.surface_rates[p] = Scalar{0.0};
1102 }
1103 }
1104 if (old_well_operable) {
1105 const std::string ctx = iterCtx.inLocalSolve() ? " (NLDD domain solve)" : "";
1106 deferred_logger.debug(" well " + this->name() + " gets STOPPED during iteration" + ctx);
1107 changed_to_stopped_this_step_ = true;
1108 }
1109 } else if (well_state.isOpen(this->name())) {
1110 this->openWell();
1111 if (!old_well_operable) {
1112 const std::string ctx = iterCtx.inLocalSolve() ? " (NLDD domain solve)" : "";
1113 deferred_logger.debug(" well " + this->name() + " gets REVIVED during iteration" + ctx);
1114 this->changed_to_open_this_step_ = true;
1115 }
1116 }
1117 }
1118
1119 template<typename TypeTag>
1120 void
1121 WellInterface<TypeTag>::addCellRates(std::map<int, RateVector>& cellRates_) const
1122 {
1123 if(!this->operability_status_.solvable)
1124 return;
1125
1126 for (int perfIdx = 0; perfIdx < this->number_of_local_perforations_; ++perfIdx) {
1127 const auto cellIdx = this->cells()[perfIdx];
1128 const auto it = cellRates_.find(cellIdx);
1129 RateVector rates = (it == cellRates_.end()) ? 0.0 : it->second;
1130 for (auto i=0*RateVector::dimension; i < RateVector::dimension; ++i)
1131 {
1132 rates[i] += connectionRates_[perfIdx][i];
1133 }
1134 cellRates_.insert_or_assign(cellIdx, rates);
1135 }
1136 }
1137
1138 template<typename TypeTag>
1141 {
1142 for (int perfIdx = 0; perfIdx < this->number_of_local_perforations_; ++perfIdx) {
1143 if (this->cells()[perfIdx] == cellIdx) {
1144 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1145 return connectionRates_[perfIdx][activeCompIdx].value();
1146 }
1147 }
1148 // this is not thread safe
1149 OPM_THROW(std::invalid_argument, "The well with name " + this->name()
1150 + " does not perforate cell " + std::to_string(cellIdx));
1151 return 0.0;
1152 }
1153
1154
1155
1156
1157 template<typename TypeTag>
1158 void
1160 checkWellOperability(const Simulator& simulator,
1161 const WellStateType& well_state,
1162 const GroupStateHelperType& groupStateHelper)
1163 {
1164 auto& deferred_logger = groupStateHelper.deferredLogger();
1165 OPM_TIMEFUNCTION();
1166 if (!this->param_.check_well_operability_) {
1167 return;
1168 }
1169
1170 if (this->wellIsStopped() && !changed_to_stopped_this_step_) {
1171 return;
1172 }
1173
1174 updateWellOperability(simulator, well_state, groupStateHelper);
1175 if (!this->operability_status_.isOperableAndSolvable()) {
1176 this->operability_status_.use_vfpexplicit = true;
1177 deferred_logger.debug("EXPLICIT_LOOKUP_VFP",
1178 "well not operable, trying with explicit vfp lookup: " + this->name());
1179 updateWellOperability(simulator, well_state, groupStateHelper);
1180 }
1181 }
1182
1183
1184
1185 template<typename TypeTag>
1186 void
1189 WellStateType& well_state,
1190 const GroupState<Scalar>& group_state,
1191 GLiftEclWells& ecl_well_map,
1192 DeferredLogger& deferred_logger)
1193 {
1194 OPM_TIMEFUNCTION();
1195 const auto& summary_state = simulator.vanguard().summaryState();
1196 const auto& well_name = this->name();
1197 if (!this->wellHasTHPConstraints(summary_state)) {
1198 const std::string msg = fmt::format("GLIFT WTEST: Well {} does not have THP constraints", well_name);
1199 deferred_logger.info(msg);
1200 return;
1201 }
1202 const auto& schedule = simulator.vanguard().schedule();
1203 const auto report_step_idx = simulator.episodeIndex();
1204 const auto& glo = schedule.glo(report_step_idx);
1205 if (!glo.has_well(well_name)) {
1206 const std::string msg = fmt::format(
1207 "GLIFT WTEST: Well {} : Gas lift not activated: "
1208 "WLIFTOPT is probably missing. Skipping.", well_name);
1209 deferred_logger.info(msg);
1210 return;
1211 }
1212 const auto& gl_well = glo.well(well_name);
1213
1214 // Use gas lift optimization to get ALQ for well test
1215 std::unique_ptr<GasLiftSingleWell> glift =
1216 initializeGliftWellTest_<GasLiftSingleWell>(simulator,
1217 well_state,
1218 group_state,
1219 ecl_well_map,
1220 deferred_logger);
1221 auto [wtest_alq, success] = glift->wellTestALQ();
1222 std::string msg;
1223 const auto& unit_system = schedule.getUnits();
1224 if (success) {
1225 well_state.well(well_name).alq_state.set(wtest_alq);
1226 msg = fmt::format(
1227 "GLIFT WTEST: Well {} : Setting ALQ to optimized value = {}",
1228 well_name, unit_system.from_si(UnitSystem::measure::gas_surface_rate, wtest_alq));
1229 }
1230 else {
1231 if (!gl_well.use_glo()) {
1232 msg = fmt::format(
1233 "GLIFT WTEST: Well {} : Gas lift optimization deactivated. Setting ALQ to WLIFTOPT item 3 = {}",
1234 well_name,
1235 unit_system.from_si(UnitSystem::measure::gas_surface_rate, well_state.well(well_name).alq_state.get()));
1236
1237 }
1238 else {
1239 msg = fmt::format(
1240 "GLIFT WTEST: Well {} : Gas lift optimization failed, no ALQ set.",
1241 well_name);
1242 }
1243 }
1244 deferred_logger.info(msg);
1245 }
1246
1247 template<typename TypeTag>
1248 void
1250 updateWellOperability(const Simulator& simulator,
1251 const WellStateType& well_state,
1252 const GroupStateHelperType& groupStateHelper)
1253 {
1254 auto& deferred_logger = groupStateHelper.deferredLogger();
1255 OPM_TIMEFUNCTION();
1256 if (this->param_.local_well_solver_control_switching_) {
1257 const bool success = updateWellOperabilityFromWellEq(simulator, groupStateHelper);
1258 if (!success) {
1259 this->operability_status_.solvable = false;
1260 deferred_logger.debug("Operability check using well equations did not converge for well "
1261 + this->name() + ". Mark the well as unsolvable." );
1262 }
1263 return;
1264 }
1265 this->operability_status_.resetOperability();
1266
1267 bool thp_controlled = this->isInjector() ? well_state.well(this->index_of_well_).injection_cmode == Well::InjectorCMode::THP:
1268 well_state.well(this->index_of_well_).production_cmode == Well::ProducerCMode::THP;
1269 bool bhp_controlled = this->isInjector() ? well_state.well(this->index_of_well_).injection_cmode == Well::InjectorCMode::BHP:
1270 well_state.well(this->index_of_well_).production_cmode == Well::ProducerCMode::BHP;
1271
1272 // Operability checking is not free
1273 // Only check wells under BHP and THP control
1274 bool check_thp = thp_controlled || this->operability_status_.thp_limit_violated_but_not_switched;
1275 if (check_thp || bhp_controlled) {
1276 updateIPR(simulator, deferred_logger);
1277 checkOperabilityUnderBHPLimit(well_state, simulator, deferred_logger);
1278 }
1279 // we do some extra checking for wells under THP control.
1280 if (check_thp) {
1281 checkOperabilityUnderTHPLimit(simulator, well_state, groupStateHelper);
1282 }
1283 }
1284
1285 template<typename TypeTag>
1286 bool
1289 const GroupStateHelperType& groupStateHelper)
1290 {
1291 OPM_TIMEFUNCTION();
1292 // only makes sense if we're using this parameter is true
1293 assert(this->param_.local_well_solver_control_switching_);
1294 this->operability_status_.resetOperability();
1295 GroupStateHelperType groupStateHelper_copy = groupStateHelper;
1296 WellStateType well_state_copy = groupStateHelper_copy.wellState();
1297 const double dt = simulator.timeStepSize();
1298 // Ensure that groupStateHelper uses well_state_copy as WellState for iterateWellEquations()
1299 // and the guard ensures that the original well state is restored at scope exit, i.e. at
1300 // the end of this function.
1301 auto guard = groupStateHelper_copy.pushWellState(well_state_copy);
1302 // equations should be converged at this stage, so only one it is needed
1303 bool converged = iterateWellEquations(simulator, dt, groupStateHelper_copy, well_state_copy);
1304 return converged;
1305 }
1306
1307 template<typename TypeTag>
1308 void
1310 scaleSegmentRatesAndPressure([[maybe_unused]] WellStateType& well_state) const
1311 {
1312 // only relevant for MSW
1313 }
1314
1315 template<typename TypeTag>
1316 void
1318 updateWellStateWithTarget(const Simulator& simulator,
1319 const GroupStateHelperType& groupStateHelper,
1320 WellStateType& well_state) const
1321 {
1322 OPM_TIMEFUNCTION();
1323 auto& deferred_logger = groupStateHelper.deferredLogger();
1324 // only bhp and wellRates are used to initilize the primaryvariables for standard wells
1325 const auto& well = this->well_ecl_;
1326 const int well_index = this->index_of_well_;
1327 auto& ws = well_state.well(well_index);
1328 const int np = well_state.numPhases();
1329 const auto& summaryState = simulator.vanguard().summaryState();
1330 const auto& schedule = simulator.vanguard().schedule();
1331
1332 // Discard old primary variables, the new well state
1333 // may not be anywhere near the old one.
1334 ws.primaryvar.resize(0);
1335
1336 if (this->wellIsStopped()) {
1337 for (int p = 0; p<np; ++p) {
1338 ws.surface_rates[p] = 0;
1339 }
1340 ws.thp = 0;
1341 return;
1342 }
1343
1344 if (this->isInjector() )
1345 {
1346 const auto& controls = well.injectionControls(summaryState);
1347
1348 InjectorType injectorType = controls.injector_type;
1349 int phasePos;
1350 switch (injectorType) {
1351 case InjectorType::WATER:
1352 {
1353 phasePos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
1354 break;
1355 }
1356 case InjectorType::OIL:
1357 {
1358 phasePos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
1359 break;
1360 }
1361 case InjectorType::GAS:
1362 {
1363 phasePos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1364 break;
1365 }
1366 default:
1367 OPM_DEFLOG_THROW(std::runtime_error, "Expected WATER, OIL or GAS as type for injectors " + this->name(), deferred_logger );
1368 }
1369
1370 const auto current = ws.injection_cmode;
1371
1372 switch (current) {
1373 case Well::InjectorCMode::RATE:
1374 {
1375 ws.surface_rates[phasePos] = (1.0 - this->rsRvInj()) * controls.surface_rate;
1376 if(this->rsRvInj() > 0) {
1377 if (injectorType == InjectorType::OIL && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1378 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1379 ws.surface_rates[gas_pos] = controls.surface_rate * this->rsRvInj();
1380 } else if (injectorType == InjectorType::GAS && FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1381 const int oil_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
1382 ws.surface_rates[oil_pos] = controls.surface_rate * this->rsRvInj();
1383 } else {
1384 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 );
1385 }
1386 }
1387 break;
1388 }
1389
1390 case Well::InjectorCMode::RESV:
1391 {
1392 std::vector<Scalar> convert_coeff(this->number_of_phases_, 1.0);
1393 this->rateConverter_.calcCoeff(/*fipreg*/ 0, this->pvtRegionIdx_, convert_coeff);
1394 const Scalar coeff = convert_coeff[phasePos];
1395 ws.surface_rates[phasePos] = controls.reservoir_rate/coeff;
1396 break;
1397 }
1398
1399 case Well::InjectorCMode::THP:
1400 {
1401 auto rates = ws.surface_rates;
1402 Scalar bhp = WellBhpThpCalculator(*this).calculateBhpFromThp(well_state,
1403 rates,
1404 well,
1405 summaryState,
1406 this->getRefDensity(),
1407 deferred_logger);
1408 ws.bhp = bhp;
1409 ws.thp = this->getTHPConstraint(summaryState);
1410
1411 // if the total rates are negative or zero
1412 // we try to provide a better intial well rate
1413 // using the well potentials
1414 Scalar total_rate = std::accumulate(rates.begin(), rates.end(), 0.0);
1415 if (total_rate <= 0.0)
1416 ws.surface_rates = ws.well_potentials;
1417
1418 break;
1419 }
1420 case Well::InjectorCMode::BHP:
1421 {
1422 ws.bhp = controls.bhp_limit;
1423 Scalar total_rate = 0.0;
1424 for (int p = 0; p<np; ++p) {
1425 total_rate += ws.surface_rates[p];
1426 }
1427 // if the total rates are negative or zero
1428 // we try to provide a better intial well rate
1429 // using the well potentials
1430 if (total_rate <= 0.0)
1431 ws.surface_rates = ws.well_potentials;
1432
1433 break;
1434 }
1435 case Well::InjectorCMode::GRUP:
1436 {
1437 assert(well.isAvailableForGroupControl());
1438 const auto& group = schedule.getGroup(well.groupName(), this->currentStep());
1439 const Scalar efficiencyFactor = well.getEfficiencyFactor() *
1440 well_state[well.name()].efficiency_scaling_factor;
1441 std::optional<Scalar> target =
1442 this->getGroupInjectionTargetRate(group,
1443 groupStateHelper,
1444 injectorType,
1445 efficiencyFactor);
1446 if (target)
1447 ws.surface_rates[phasePos] = *target;
1448 break;
1449 }
1450 case Well::InjectorCMode::CMODE_UNDEFINED:
1451 {
1452 OPM_DEFLOG_THROW(std::runtime_error, "Well control must be specified for well " + this->name(), deferred_logger );
1453 }
1454
1455 }
1456 // for wells with zero injection rate, if we assign exactly zero rate,
1457 // we will have to assume some trivial composition in the wellbore.
1458 // here, we use some small value (about 0.01 m^3/day ~= 1.e-7) to initialize
1459 // the zero rate target, then we can use to retain the composition information
1460 // within the wellbore from the previous result, and hopefully it is a good
1461 // initial guess for the zero rate target.
1462 ws.surface_rates[phasePos] = std::max(Scalar{1.e-7}, ws.surface_rates[phasePos]);
1463
1464 if (ws.bhp == 0.) {
1465 ws.bhp = controls.bhp_limit;
1466 }
1467 }
1468 //Producer
1469 else
1470 {
1471 const auto current = ws.production_cmode;
1472 const auto& controls = well.productionControls(summaryState);
1473 switch (current) {
1474 case Well::ProducerCMode::ORAT:
1475 {
1476 const int oil_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
1477 Scalar current_rate = -ws.surface_rates[oil_pos];
1478 // for trivial rates or opposite direction we don't just scale the rates
1479 // but use either the potentials or the mobility ratio to initial the well rates
1480 if (current_rate > 0.0) {
1481 for (int p = 0; p<np; ++p) {
1482 ws.surface_rates[p] *= controls.oil_rate/current_rate;
1483 }
1484 } else {
1485 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1486 double control_fraction = fractions[oil_pos];
1487 if (control_fraction != 0.0) {
1488 for (int p = 0; p<np; ++p) {
1489 ws.surface_rates[p] = - fractions[p] * controls.oil_rate/control_fraction;
1490 }
1491 }
1492 }
1493 break;
1494 }
1495 case Well::ProducerCMode::WRAT:
1496 {
1497 const int water_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
1498 Scalar current_rate = -ws.surface_rates[water_pos];
1499 // for trivial rates or opposite direction we don't just scale the rates
1500 // but use either the potentials or the mobility ratio to initial the well rates
1501 if (current_rate > 0.0) {
1502 for (int p = 0; p<np; ++p) {
1503 ws.surface_rates[p] *= controls.water_rate/current_rate;
1504 }
1505 } else {
1506 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1507 const Scalar control_fraction = fractions[water_pos];
1508 if (control_fraction != 0.0) {
1509 for (int p = 0; p<np; ++p) {
1510 ws.surface_rates[p] = - fractions[p] * controls.water_rate / control_fraction;
1511 }
1512 }
1513 }
1514 break;
1515 }
1516 case Well::ProducerCMode::GRAT:
1517 {
1518 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1519 Scalar current_rate = -ws.surface_rates[gas_pos];
1520 // or trivial rates or opposite direction we don't just scale the rates
1521 // but use either the potentials or the mobility ratio to initial the well rates
1522 if (current_rate > 0.0) {
1523 for (int p = 0; p<np; ++p) {
1524 ws.surface_rates[p] *= controls.gas_rate/current_rate;
1525 }
1526 } else {
1527 const std::vector<Scalar > fractions = initialWellRateFractions(simulator, well_state);
1528 const Scalar control_fraction = fractions[gas_pos];
1529 if (control_fraction != 0.0) {
1530 for (int p = 0; p<np; ++p) {
1531 ws.surface_rates[p] = - fractions[p] * controls.gas_rate / control_fraction;
1532 }
1533 }
1534 }
1535
1536 break;
1537
1538 }
1539 case Well::ProducerCMode::LRAT:
1540 {
1541 const int water_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
1542 const int oil_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
1543 Scalar current_rate = - ws.surface_rates[water_pos]
1544 - ws.surface_rates[oil_pos];
1545 // or trivial rates or opposite direction we don't just scale the rates
1546 // but use either the potentials or the mobility ratio to initial the well rates
1547 if (current_rate > 0.0) {
1548 for (int p = 0; p<np; ++p) {
1549 ws.surface_rates[p] *= controls.liquid_rate/current_rate;
1550 }
1551 } else {
1552 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1553 const Scalar control_fraction = fractions[water_pos] + fractions[oil_pos];
1554 if (control_fraction != 0.0) {
1555 for (int p = 0; p<np; ++p) {
1556 ws.surface_rates[p] = - fractions[p] * controls.liquid_rate / control_fraction;
1557 }
1558 }
1559 }
1560 break;
1561 }
1562 case Well::ProducerCMode::CRAT:
1563 {
1564 OPM_DEFLOG_THROW(std::runtime_error,
1565 fmt::format("CRAT control not supported, well {}", this->name()),
1566 deferred_logger);
1567 }
1568 case Well::ProducerCMode::RESV:
1569 {
1570 std::vector<Scalar> convert_coeff(this->number_of_phases_, 1.0);
1571 this->rateConverter_.calcCoeff(/*fipreg*/ 0, this->pvtRegionIdx_, ws.surface_rates, convert_coeff);
1572 Scalar total_res_rate = 0.0;
1573 for (int p = 0; p<np; ++p) {
1574 total_res_rate -= ws.surface_rates[p] * convert_coeff[p];
1575 }
1576 if (controls.prediction_mode) {
1577 // or trivial rates or opposite direction we don't just scale the rates
1578 // but use either the potentials or the mobility ratio to initial the well rates
1579 if (total_res_rate > 0.0) {
1580 for (int p = 0; p<np; ++p) {
1581 ws.surface_rates[p] *= controls.resv_rate/total_res_rate;
1582 }
1583 } else {
1584 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1585 for (int p = 0; p<np; ++p) {
1586 ws.surface_rates[p] = - fractions[p] * controls.resv_rate / convert_coeff[p];
1587 }
1588 }
1589 } else {
1590 std::vector<Scalar> hrates(this->number_of_phases_,0.);
1591 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
1592 const int phase_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
1593 hrates[phase_pos] = controls.water_rate;
1594 }
1595 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1596 const int phase_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
1597 hrates[phase_pos] = controls.oil_rate;
1598 }
1599 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1600 const int phase_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
1601 hrates[phase_pos] = controls.gas_rate;
1602 }
1603 std::vector<Scalar> hrates_resv(this->number_of_phases_,0.);
1604 this->rateConverter_.calcReservoirVoidageRates(/*fipreg*/ 0, this->pvtRegionIdx_, hrates, hrates_resv);
1605 Scalar target = std::accumulate(hrates_resv.begin(), hrates_resv.end(), 0.0);
1606 // or trivial rates or opposite direction we don't just scale the rates
1607 // but use either the potentials or the mobility ratio to initial the well rates
1608 if (total_res_rate > 0.0) {
1609 for (int p = 0; p<np; ++p) {
1610 ws.surface_rates[p] *= target/total_res_rate;
1611 }
1612 } else {
1613 const std::vector<Scalar> fractions = initialWellRateFractions(simulator, well_state);
1614 for (int p = 0; p<np; ++p) {
1615 ws.surface_rates[p] = - fractions[p] * target / convert_coeff[p];
1616 }
1617 }
1618 }
1619 break;
1620 }
1621 case Well::ProducerCMode::BHP:
1622 {
1623 ws.bhp = controls.bhp_limit;
1624 Scalar total_rate = 0.0;
1625 for (int p = 0; p<np; ++p) {
1626 total_rate -= ws.surface_rates[p];
1627 }
1628 // if the total rates are negative or zero
1629 // we try to provide a better intial well rate
1630 // using the well potentials
1631 if (total_rate <= 0.0){
1632 for (int p = 0; p<np; ++p) {
1633 ws.surface_rates[p] = -ws.well_potentials[p];
1634 }
1635 }
1636 break;
1637 }
1638 case Well::ProducerCMode::THP:
1639 {
1640 const bool update_success = updateWellStateWithTHPTargetProd(simulator, well_state, groupStateHelper);
1641
1642 if (!update_success) {
1643 // the following is the original way of initializing well state with THP constraint
1644 // keeping it for robust reason in case that it fails to get a bhp value with THP constraint
1645 // more sophisticated design might be needed in the future
1646 auto rates = ws.surface_rates;
1647 this->adaptRatesForVFP(rates);
1649 well_state, rates, well, summaryState, this->getRefDensity(), deferred_logger);
1650 ws.bhp = bhp;
1651 ws.thp = this->getTHPConstraint(summaryState);
1652 // if the total rates are negative or zero
1653 // we try to provide a better initial well rate
1654 // using the well potentials
1655 const Scalar total_rate = -std::accumulate(rates.begin(), rates.end(), 0.0);
1656 if (total_rate <= 0.0) {
1657 for (int p = 0; p < this->number_of_phases_; ++p) {
1658 ws.surface_rates[p] = -ws.well_potentials[p];
1659 }
1660 }
1661 }
1662 break;
1663 }
1664 case Well::ProducerCMode::GRUP:
1665 {
1666 assert(well.isAvailableForGroupControl());
1667 const auto& group = schedule.getGroup(well.groupName(), this->currentStep());
1668 const Scalar efficiencyFactor = well.getEfficiencyFactor() *
1669 well_state[well.name()].efficiency_scaling_factor;
1670 Scalar scale = this->getGroupProductionTargetRate(group,
1671 groupStateHelper,
1672 efficiencyFactor);
1673
1674 // we don't want to scale with zero and get zero rates.
1675 if (scale > 0) {
1676 for (int p = 0; p<np; ++p) {
1677 ws.surface_rates[p] *= scale;
1678 }
1679 ws.trivial_group_target = false;
1680 } else {
1681 // If group target is trivial we dont want to flip to other controls. To avoid oscillation we store
1682 // this information in the well state and explicitly check for this condition when evaluating well controls.
1683 ws.trivial_group_target = true;
1684 }
1685 break;
1686 }
1687 case Well::ProducerCMode::CMODE_UNDEFINED:
1689 {
1690 OPM_DEFLOG_THROW(std::runtime_error, "Well control must be specified for well " + this->name() , deferred_logger);
1691 break;
1692 }
1693 } // end of switch
1694
1695 if (ws.bhp == 0.) {
1696 ws.bhp = controls.bhp_limit;
1697 }
1698 }
1699 }
1700
1701 template<typename TypeTag>
1702 bool
1704 wellUnderZeroRateTarget(const GroupStateHelperType& groupStateHelper) const
1705 {
1706 OPM_TIMEFUNCTION();
1707 const auto& well_state = groupStateHelper.wellState();
1708 // Check if well is under zero rate control, either directly or from group
1709 const bool isGroupControlled = this->wellUnderGroupControl(well_state.well(this->index_of_well_));
1710 if (!isGroupControlled) {
1711 // well is not under group control, check "individual" version
1712 const auto& summaryState = groupStateHelper.summaryState();
1713 return this->wellUnderZeroRateTargetIndividual(summaryState, well_state);
1714 } else {
1715 return this->wellUnderZeroGroupRateTarget(groupStateHelper, isGroupControlled);
1716 }
1717 }
1718
1719 template <typename TypeTag>
1720 bool
1722 const std::optional<bool> group_control) const
1723 {
1724 const auto& well_state = groupStateHelper.wellState();
1725 // Check if well is under zero rate target from group
1726 const bool isGroupControlled = group_control.value_or(this->wellUnderGroupControl(well_state.well(this->index_of_well_)));
1727 if (isGroupControlled) {
1728 return this->zeroGroupRateTarget(groupStateHelper);
1729 }
1730 return false;
1731 }
1732
1733 template<typename TypeTag>
1734 bool
1736 stoppedOrZeroRateTarget(const GroupStateHelperType& groupStateHelper) const
1737 {
1738 // Check if well is stopped or under zero rate control, either
1739 // directly or from group.
1740 return this->wellIsStopped()
1741 || this->wellUnderZeroRateTarget(groupStateHelper);
1742 }
1743
1744 template<typename TypeTag>
1745 std::vector<typename WellInterface<TypeTag>::Scalar>
1747 initialWellRateFractions(const Simulator& simulator,
1748 const WellStateType& well_state) const
1749 {
1750 OPM_TIMEFUNCTION();
1751 const int np = this->number_of_phases_;
1752 std::vector<Scalar> scaling_factor(np);
1753 const auto& ws = well_state.well(this->index_of_well_);
1754
1755 Scalar total_potentials = 0.0;
1756 for (int p = 0; p<np; ++p) {
1757 total_potentials += ws.well_potentials[p];
1758 }
1759 if (total_potentials > 0) {
1760 for (int p = 0; p<np; ++p) {
1761 scaling_factor[p] = ws.well_potentials[p] / total_potentials;
1762 }
1763 return scaling_factor;
1764 }
1765 // if we don't have any potentials we weight it using the mobilites
1766 // We only need approximation so we don't bother with the vapporized oil and dissolved gas
1767 Scalar total_tw = 0;
1768 const int nperf = this->number_of_local_perforations_;
1769 for (int perf = 0; perf < nperf; ++perf) {
1770 total_tw += this->well_index_[perf];
1771 }
1772 total_tw = this->parallelWellInfo().communication().sum(total_tw);
1773
1774 for (int perf = 0; perf < nperf; ++perf) {
1775 const int cell_idx = this->well_cells_[perf];
1776 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
1777 const auto& fs = intQuants.fluidState();
1778 const Scalar well_tw_fraction = this->well_index_[perf] / total_tw;
1779 Scalar total_mobility = 0.0;
1780 for (int p = 0; p < np; ++p) {
1781 const int canonical_phase_idx = FluidSystem::activeToCanonicalPhaseIdx(p);
1782 total_mobility += fs.invB(canonical_phase_idx).value() * intQuants.mobility(canonical_phase_idx).value();
1783 }
1784 for (int p = 0; p < np; ++p) {
1785 const int canonical_phase_idx = FluidSystem::activeToCanonicalPhaseIdx(p);
1786 scaling_factor[p] += well_tw_fraction * fs.invB(canonical_phase_idx).value() * intQuants.mobility(canonical_phase_idx).value() / total_mobility;
1787 }
1788 }
1789 return scaling_factor;
1790 }
1791
1792
1793
1794 template <typename TypeTag>
1795 void
1798 WellStateType& well_state,
1799 DeferredLogger& deferred_logger) const
1800 {
1801 assert(this->isProducer());
1802 OPM_TIMEFUNCTION();
1803 // Check if the rates of this well only are single-phase, do nothing
1804 // if more than one nonzero rate.
1805 auto& ws = well_state.well(this->index_of_well_);
1806 int nonzero_rate_index = -1;
1807 const Scalar floating_point_error_epsilon = 1e-14;
1808 for (int p = 0; p < this->number_of_phases_; ++p) {
1809 if (std::abs(ws.surface_rates[p]) > floating_point_error_epsilon) {
1810 if (nonzero_rate_index == -1) {
1811 nonzero_rate_index = p;
1812 } else {
1813 // More than one nonzero rate.
1814 return;
1815 }
1816 }
1817 }
1818
1819 // Calculate rates at bhp limit, or 1 bar if no limit.
1820 std::vector<Scalar> well_q_s(this->number_of_phases_, 0.0);
1821 bool rates_evaluated_at_1bar = false;
1822 {
1823 const auto& summary_state = simulator.vanguard().summaryState();
1824 const auto& prod_controls = this->well_ecl_.productionControls(summary_state);
1825 const double bhp_limit = std::max(prod_controls.bhp_limit, 1.0 * unit::barsa);
1826 this->computeWellRatesWithBhp(simulator, bhp_limit, well_q_s, deferred_logger);
1827 // Remember of we evaluated the rates at (approx.) 1 bar or not.
1828 rates_evaluated_at_1bar = (bhp_limit < 1.1 * unit::barsa);
1829 // Check that no rates are positive.
1830 if (std::ranges::any_of(well_q_s, [](Scalar q) { return q > 0.0; })) {
1831 // Did we evaluate at 1 bar? If not, then we can try again at 1 bar.
1832 if (!rates_evaluated_at_1bar) {
1833 this->computeWellRatesWithBhp(simulator, 1.0 * unit::barsa, well_q_s, deferred_logger);
1834 rates_evaluated_at_1bar = true;
1835 }
1836 // At this point we can only set the wrong-direction (if any) values to zero.
1837 for (auto& q : well_q_s) {
1838 q = std::min(q, Scalar{0.0});
1839 }
1840 }
1841 }
1842
1843 if (nonzero_rate_index == -1) {
1844 // No nonzero rates on input.
1845 // Use the computed rate directly, or scaled by a factor
1846 // 0.5 (to avoid too high values) if it was evaluated at 1 bar.
1847 const Scalar factor = rates_evaluated_at_1bar ? 0.5 : 1.0;
1848 for (int p = 0; p < this->number_of_phases_; ++p) {
1849 ws.surface_rates[p] = factor * well_q_s[p];
1850 }
1851 return;
1852 }
1853
1854 // If we are here, we had a single nonzero rate for the well,
1855 // typically from a rate constraint. We must make sure it is
1856 // respected, so if it was lower than the calculated rate for
1857 // the same phase we scale all rates to match.
1858 const Scalar initial_nonzero_rate = ws.surface_rates[nonzero_rate_index];
1859 const Scalar computed_rate = well_q_s[nonzero_rate_index];
1860 if (std::abs(initial_nonzero_rate) < std::abs(computed_rate)) {
1861 // Note that both rates below are negative. The factor should be < 1.0.
1862 const Scalar factor = initial_nonzero_rate / computed_rate;
1863 assert(factor < 1.0);
1864 for (int p = 0; p < this->number_of_phases_; ++p) {
1865 // We skip the nonzero_rate_index, as that should remain as it was.
1866 if (p != nonzero_rate_index) {
1867 ws.surface_rates[p] = factor * well_q_s[p];
1868 }
1869 }
1870 return;
1871 }
1872
1873 // If we are here, we had a single nonzero rate, but it was
1874 // higher than the one calculated from the bhp limit, so we
1875 // use the calculated rates.
1876 for (int p = 0; p < this->number_of_phases_; ++p) {
1877 ws.surface_rates[p] = well_q_s[p];
1878 }
1879 }
1880
1881 template <typename TypeTag>
1882 template<class Value>
1883 void
1885 getTw(std::vector<Value>& Tw,
1886 const int perf,
1887 const IntensiveQuantities& intQuants,
1888 const Value& trans_mult,
1889 const SingleWellStateType& ws) const
1890 {
1891 OPM_TIMEFUNCTION_LOCAL(Subsystem::Wells);
1892 // Add a Forchheimer term to the gas phase CTF if the run uses
1893 // either of the WDFAC or the WDFACCOR keywords.
1894 if (static_cast<std::size_t>(perf) >= this->well_cells_.size()) {
1895 OPM_THROW(std::invalid_argument,"The perforation index exceeds the size of the local containers - possibly wellIndex was called with a global instead of a local perforation index!");
1896 }
1897
1898 if constexpr (! Indices::gasEnabled) {
1899 return;
1900 }
1901
1902 const auto& wdfac = this->well_ecl_.getWDFAC();
1903
1904 if (! wdfac.useDFactor() || (this->well_index_[perf] == 0.0)) {
1905 return;
1906 }
1907
1908 const Scalar d = this->computeConnectionDFactor(perf, intQuants, ws);
1909 if (d < 1.0e-15) {
1910 return;
1911 }
1912
1913 // Solve quadratic equations for connection rates satisfying the ipr and the flow-dependent skin.
1914 // If more than one solution, pick the one corresponding to lowest absolute rate (smallest skin).
1915 const auto& connection = this->well_ecl_.getConnections()[ws.perf_data.ecl_index[perf]];
1916 const Scalar Kh = connection.Kh();
1917 const Scalar scaling = std::numbers::pi * Kh * connection.wpimult();
1918 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1919
1920 const Scalar connection_pressure = ws.perf_data.pressure[perf];
1921 const Scalar cell_pressure = getValue(intQuants.fluidState().pressure(FluidSystem::gasPhaseIdx));
1922 const Scalar drawdown = cell_pressure - connection_pressure;
1923 const Scalar invB = getValue(intQuants.fluidState().invB(FluidSystem::gasPhaseIdx));
1924 const Scalar mob_g = getValue(intQuants.mobility(FluidSystem::gasPhaseIdx)) * invB;
1925 const Scalar a = d;
1926 const Scalar b = 2 * scaling / getValue(Tw[gas_comp_idx]);
1927 const Scalar c = -2 * scaling * mob_g * drawdown;
1928
1929 Scalar consistent_Q = -1.0e20;
1930 // Find and check negative solutions (a --> -a)
1931 const Scalar r2n = b*b + 4*a*c;
1932 if (r2n >= 0) {
1933 const Scalar rn = std::sqrt(r2n);
1934 const Scalar xn1 = (b-rn)*0.5/a;
1935 if (xn1 <= 0) {
1936 consistent_Q = xn1;
1937 }
1938 const Scalar xn2 = (b+rn)*0.5/a;
1939 if (xn2 <= 0 && xn2 > consistent_Q) {
1940 consistent_Q = xn2;
1941 }
1942 }
1943 // Find and check positive solutions
1944 consistent_Q *= -1;
1945 const Scalar r2p = b*b - 4*a*c;
1946 if (r2p >= 0) {
1947 const Scalar rp = std::sqrt(r2p);
1948 const Scalar xp1 = (rp-b)*0.5/a;
1949 if (xp1 > 0 && xp1 < consistent_Q) {
1950 consistent_Q = xp1;
1951 }
1952 const Scalar xp2 = -(rp+b)*0.5/a;
1953 if (xp2 > 0 && xp2 < consistent_Q) {
1954 consistent_Q = xp2;
1955 }
1956 }
1957 Tw[gas_comp_idx] = 1.0 / (1.0 / (trans_mult * this->well_index_[perf]) + (consistent_Q/2 * d / scaling));
1958 }
1959
1960 template <typename TypeTag>
1961 void
1963 updateConnectionDFactor(const Simulator& simulator,
1964 SingleWellStateType& ws) const
1965 {
1966 if (! this->well_ecl_.getWDFAC().useDFactor()) {
1967 return;
1968 }
1969
1970 auto& d_factor = ws.perf_data.connection_d_factor;
1971
1972 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
1973 const int cell_idx = this->well_cells_[perf];
1974 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1975
1976 d_factor[perf] = this->computeConnectionDFactor(perf, intQuants, ws);
1977 }
1978 }
1979
1980 template <typename TypeTag>
1983 computeConnectionDFactor(const int perf,
1984 const IntensiveQuantities& intQuants,
1985 const SingleWellStateType& ws) const
1986 {
1987 auto rhoGS = [regIdx = this->pvtRegionIdx()]() {
1988 return FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, regIdx);
1989 };
1990
1991 // Viscosity is evaluated at connection pressure.
1992 auto gas_visc = [connection_pressure = ws.perf_data.pressure[perf],
1993 temperature = ws.temperature,
1994 regIdx = this->pvtRegionIdx(), &intQuants]()
1995 {
1996 const auto rv = getValue(intQuants.fluidState().Rv());
1997
1998 const auto& gasPvt = FluidSystem::gasPvt();
1999
2000 // Note that rv here is from grid block with typically
2001 // p_block > connection_pressure
2002 // so we may very well have rv > rv_sat
2003 const Scalar rv_sat = gasPvt.saturatedOilVaporizationFactor
2004 (regIdx, temperature, connection_pressure);
2005
2006 if (! (rv < rv_sat)) {
2007 return gasPvt.saturatedViscosity(regIdx, temperature,
2008 connection_pressure);
2009 }
2010
2011 return gasPvt.viscosity(regIdx, temperature, connection_pressure,
2012 rv, getValue(intQuants.fluidState().Rvw()));
2013 };
2014
2015 const auto& connection = this->well_ecl_.getConnections()
2016 [ws.perf_data.ecl_index[perf]];
2017
2018 return this->well_ecl_.getWDFAC().getDFactor(rhoGS, gas_visc, connection);
2019 }
2020
2021
2022 template <typename TypeTag>
2023 void
2026 SingleWellStateType& ws) const
2027 {
2028 auto connCF = [&connIx = std::as_const(ws.perf_data.ecl_index),
2029 &conns = this->well_ecl_.getConnections()]
2030 (const int perf)
2031 {
2032 return conns[connIx[perf]].CF();
2033 };
2034
2035 auto obtain = [](const Eval& value)
2036 {
2037 return getValue(value);
2038 };
2039
2040 auto& tmult = ws.perf_data.connection_compaction_tmult;
2041 auto& ctf = ws.perf_data.connection_transmissibility_factor;
2042
2043 for (int perf = 0; perf < this->number_of_local_perforations_; ++perf) {
2044 const int cell_idx = this->well_cells_[perf];
2045 Scalar trans_mult(0.0);
2046 getTransMult(trans_mult, simulator, cell_idx, obtain);
2047 tmult[perf] = trans_mult;
2048
2049 ctf[perf] = connCF(perf) * tmult[perf];
2050 }
2051 }
2052
2053
2054 template<typename TypeTag>
2057 {
2058 if constexpr (Indices::oilEnabled) {
2059 return fs.pressure(FluidSystem::oilPhaseIdx);
2060 } else if constexpr (Indices::gasEnabled) {
2061 return fs.pressure(FluidSystem::gasPhaseIdx);
2062 } else {
2063 return fs.pressure(FluidSystem::waterPhaseIdx);
2064 }
2065 }
2066
2067 template <typename TypeTag>
2068 template<class Value, class Callback>
2069 void
2071 getTransMult(Value& trans_mult,
2072 const Simulator& simulator,
2073 const int cell_idx,
2074 Callback& extendEval) const
2075 {
2076 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2077 trans_mult = simulator.problem().template wellTransMultiplier<Value>(intQuants, cell_idx, extendEval);
2078 }
2079
2080 template <typename TypeTag>
2081 template<class Value, class Callback>
2082 void
2084 getMobility(const Simulator& simulator,
2085 const int local_perf_index,
2086 std::vector<Value>& mob,
2087 Callback& extendEval,
2088 [[maybe_unused]] DeferredLogger& deferred_logger) const
2089 {
2090 auto relpermArray = []()
2091 {
2092 if constexpr (std::is_same_v<Value, Scalar>) {
2093 return std::array<Scalar,3>{};
2094 } else {
2095 return std::array<Eval,3>{};
2096 }
2097 };
2098 if (static_cast<std::size_t>(local_perf_index) >= this->well_cells_.size()) {
2099 OPM_THROW(std::invalid_argument,"The perforation index exceeds the size of the local containers - possibly getMobility was called with a global instead of a local perforation index!");
2100 }
2101 const int cell_idx = this->well_cells_[local_perf_index];
2102 assert (int(mob.size()) == this->num_conservation_quantities_);
2103 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
2104 const auto& materialLawManager = simulator.problem().materialLawManager();
2105
2106 // either use mobility of the perforation cell or calculate its own
2107 // based on passing the saturation table index
2108 const int satid = this->saturation_table_number_[local_perf_index] - 1;
2109 const int satid_elem = materialLawManager->satnumRegionIdx(cell_idx);
2110 if (satid == satid_elem) { // the same saturation number is used. i.e. just use the mobilty from the cell
2111 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
2112 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2113 continue;
2114 }
2115
2116 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2117 mob[activeCompIdx] = extendEval(intQuants.mobility(phaseIdx));
2118 }
2119 if constexpr (has_solvent) {
2120 mob[Indices::contiSolventEqIdx] = extendEval(intQuants.solventMobility());
2121 }
2122 } else {
2123 const auto& paramsCell = materialLawManager->connectionMaterialLawParams(satid, cell_idx);
2124 auto relativePerms = relpermArray();
2125 MaterialLaw::relativePermeabilities(relativePerms, paramsCell, intQuants.fluidState());
2126
2127 // reset the satnumvalue back to original
2128 materialLawManager->connectionMaterialLawParams(satid_elem, cell_idx);
2129
2130 // compute the mobility
2131 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
2132 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2133 continue;
2134 }
2135
2136 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
2137 mob[activeCompIdx] = extendEval(relativePerms[phaseIdx] / intQuants.fluidState().viscosity(phaseIdx));
2138 }
2139
2140 if constexpr (has_solvent) {
2141 const auto Fsolgas = intQuants.solventSaturation() / (intQuants.solventSaturation() + intQuants.fluidState().saturation(FluidSystem::gasPhaseIdx));
2142 using SolventModule = BlackOilSolventModule<TypeTag>;
2143 if (Fsolgas > SolventModule::cutOff) { // same cutoff as in the solvent model to avoid division by zero
2144 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(FluidSystem::gasPhaseIdx));
2145 const auto& ssfnKrg = SolventModule::ssfnKrg(satid);
2146 const auto& ssfnKrs = SolventModule::ssfnKrs(satid);
2147 mob[activeGasCompIdx] *= extendEval(ssfnKrg.eval(1-Fsolgas, /*extrapolate=*/true));
2148 mob[Indices::contiSolventEqIdx] = extendEval(ssfnKrs.eval(Fsolgas, /*extrapolate=*/true) * relativePerms[activeGasCompIdx] / intQuants.solventViscosity());
2149 }
2150 }
2151 }
2152
2153 if (this->isInjector() && !this->inj_fc_multiplier_.empty()) {
2154 const auto perf_ecl_index = this->perforationData()[local_perf_index].ecl_index;
2155 const auto& connections = this->well_ecl_.getConnections();
2156 const auto& connection = connections[perf_ecl_index];
2157 if (connection.filterCakeActive()) {
2158 std::ranges::transform(mob, mob.begin(),
2159 [mult = this->inj_fc_multiplier_[local_perf_index]]
2160 (const auto val)
2161 { return val * mult; });
2162 }
2163 }
2164 }
2165
2166
2167 template<typename TypeTag>
2168 bool
2171 WellStateType& well_state,
2172 const GroupStateHelperType& groupStateHelper) const
2173 {
2174 auto& deferred_logger = groupStateHelper.deferredLogger();
2175 OPM_TIMEFUNCTION();
2176 const auto& summary_state = simulator.vanguard().summaryState();
2177
2178 auto bhp_at_thp_limit = computeBhpAtThpLimitProdWithAlq(
2179 simulator, groupStateHelper, summary_state, this->getALQ(well_state), /*iterate_if_no_solution */ false);
2180 if (bhp_at_thp_limit) {
2181 std::vector<Scalar> rates(this->number_of_phases_, 0.0);
2182 if (thp_update_iterations) {
2183 computeWellRatesWithBhpIterations(simulator, *bhp_at_thp_limit,
2184 groupStateHelper, rates);
2185 } else {
2186 computeWellRatesWithBhp(simulator, *bhp_at_thp_limit,
2187 rates, deferred_logger);
2188 }
2189 auto& ws = well_state.well(this->name());
2190 ws.surface_rates = rates;
2191 ws.bhp = *bhp_at_thp_limit;
2192 ws.thp = this->getTHPConstraint(summary_state);
2193 return true;
2194 } else {
2195 return false;
2196 }
2197 }
2198
2199 template<typename TypeTag>
2200 std::optional<typename WellInterface<TypeTag>::Scalar>
2203 const WellStateType& well_state,
2204 Scalar bhp,
2205 const SummaryState& summary_state,
2206 const Scalar alq_value)
2207 {
2208 OPM_TIMEFUNCTION();
2209 WellStateType well_state_copy = well_state;
2210 const auto& groupStateHelper = simulator.problem().wellModel().groupStateHelper();
2211 GroupStateHelperType groupStateHelper_copy = groupStateHelper;
2212 auto well_guard = groupStateHelper_copy.pushWellState(well_state_copy);
2213 const double dt = simulator.timeStepSize();
2214 const bool converged = this->solveWellWithBhp(
2215 simulator, dt, bhp, groupStateHelper_copy, well_state_copy
2216 );
2217
2218 bool zero_rates;
2219 auto rates = well_state_copy.well(this->index_of_well_).surface_rates;
2220 zero_rates = true;
2221 for (std::size_t p = 0; p < rates.size(); ++p) {
2222 zero_rates &= rates[p] == 0.0;
2223 }
2224 // For zero rates or unconverged bhp the implicit IPR is problematic.
2225 // Use the old approach for now
2226 if (zero_rates || !converged) {
2227 return this->computeBhpAtThpLimitProdWithAlq(simulator, groupStateHelper_copy, summary_state, alq_value, /*iterate_if_no_solution */ false);
2228 }
2229 this->updateIPRImplicit(simulator, groupStateHelper_copy, well_state_copy);
2230 this->adaptRatesForVFP(rates);
2231 return WellBhpThpCalculator(*this).estimateStableBhp(well_state_copy, this->well_ecl_, rates, this->getRefDensity(), summary_state, alq_value);
2232 }
2233
2234 template <typename TypeTag>
2235 void
2238 const std::function<Scalar(const Scalar)>& connPICalc,
2239 const std::vector<Scalar>& mobility,
2240 Scalar* connPI) const
2241 {
2242 const int np = this->number_of_phases_;
2243 for (int p = 0; p < np; ++p) {
2244 // Note: E100's notion of PI value phase mobility includes
2245 // the reciprocal FVF.
2246 const int canonical_phase_idx = FluidSystem::activeToCanonicalPhaseIdx(p);
2247 const auto connMob =
2248 mobility[FluidSystem::activePhaseToActiveCompIdx(p)] * fs.invB(canonical_phase_idx).value();
2249
2250 connPI[p] = connPICalc(connMob);
2251 }
2252
2253 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) &&
2254 FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx))
2255 {
2256 const auto io = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
2257 const auto ig = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
2258
2259 const auto vapoil = connPI[ig] * fs.Rv().value();
2260 const auto disgas = connPI[io] * fs.Rs().value();
2261
2262 connPI[io] += vapoil;
2263 connPI[ig] += disgas;
2264 }
2265 }
2266
2267
2268 template <typename TypeTag>
2269 void
2272 const Phase preferred_phase,
2273 const std::function<Scalar(const Scalar)>& connIICalc,
2274 const std::vector<Scalar>& mobility,
2275 Scalar* connII,
2276 DeferredLogger& deferred_logger) const
2277 {
2278 auto phase_pos = 0;
2279 if (preferred_phase == Phase::GAS) {
2280 phase_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
2281 }
2282 else if (preferred_phase == Phase::OIL) {
2283 phase_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
2284 }
2285 else if (preferred_phase == Phase::WATER) {
2286 phase_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
2287 }
2288 else {
2289 OPM_DEFLOG_THROW(NotImplemented,
2290 fmt::format("Unsupported Injector Type ({}) "
2291 "for well {} during connection I.I. calculation",
2292 static_cast<int>(preferred_phase), this->name()),
2293 deferred_logger);
2294 }
2295
2296 const auto mt = std::accumulate(mobility.begin(), mobility.end(), 0.0);
2297 const int canonicalPhaseIdx = FluidSystem::activeToCanonicalPhaseIdx(phase_pos);
2298 connII[phase_pos] = connIICalc(mt * fs.invB(canonicalPhaseIdx).value());
2299 }
2300
2301 template<typename TypeTag>
2302 template<class GasLiftSingleWell>
2303 std::unique_ptr<GasLiftSingleWell>
2305 initializeGliftWellTest_(const Simulator& simulator,
2306 WellStateType& well_state,
2307 const GroupState<Scalar>& group_state,
2308 GLiftEclWells& ecl_well_map,
2309 DeferredLogger& deferred_logger)
2310 {
2311 // Instantiate group info object (without initialization) since it is needed in GasLiftSingleWell
2312 auto& comm = simulator.vanguard().grid().comm();
2313 ecl_well_map.try_emplace(this->name(), &(this->wellEcl()), this->indexOfWell());
2314 const auto& iterCtx = simulator.problem().iterationContext();
2316 ecl_well_map,
2317 simulator.vanguard().schedule(),
2318 simulator.vanguard().summaryState(),
2319 simulator.episodeIndex(),
2320 iterCtx,
2321 deferred_logger,
2322 well_state,
2323 group_state,
2324 comm,
2325 false
2326 };
2327
2328 // Return GasLiftSingleWell object to use the wellTestALQ() function
2329 std::set<int> sync_groups;
2330 const auto& summary_state = simulator.vanguard().summaryState();
2331 return std::make_unique<GasLiftSingleWell>(*this,
2332 simulator,
2333 summary_state,
2334 deferred_logger,
2335 well_state,
2336 group_state,
2337 group_info,
2338 sync_groups,
2339 comm,
2340 false);
2341
2342 }
2343
2344} // namespace Opm
2345
2346#endif
#define OPM_DEFLOG_THROW(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:45
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:69
Definition: DeferredLogger.hpp:57
void info(const std::string &tag, const std::string &message)
Definition: GasLiftGroupInfo.hpp:47
Definition: GroupStateHelper.hpp:57
GroupState< Scalar > & groupState() const
Definition: GroupStateHelper.hpp:298
const SummaryState & summaryState() const
Definition: GroupStateHelper.hpp:415
const WellState< Scalar, IndexTraits > & wellState() const
Definition: GroupStateHelper.hpp:494
DeferredLogger & deferredLogger() const
Get the deferred logger.
Definition: GroupStateHelper.hpp:234
WellStateGuard pushWellState(WellState< Scalar, IndexTraits > &well_state)
Definition: GroupStateHelper.hpp:354
GroupStateGuard pushGroupState(GroupState< Scalar > &group_state)
Definition: GroupStateHelper.hpp:331
Definition: GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: SingleWellState.hpp:43
Scalar temperature
Definition: SingleWellState.hpp:105
PerfData< Scalar > perf_data
Definition: SingleWellState.hpp:143
Class for computing BHP limits.
Definition: WellBhpThpCalculator.hpp:41
Scalar calculateMinimumBhpFromThp(const WellState< Scalar, IndexTraits > &well_state, const Well &well, const SummaryState &summaryState, const Scalar rho) const
Scalar mostStrictBhpFromBhpLimits(const SummaryState &summaryState) const
Obtain the most strict BHP from BHP limits.
bool isStableSolution(const WellState< Scalar, IndexTraits > &well_state, const Well &well, const std::vector< Scalar > &rates, const SummaryState &summaryState) const
EvalWell calculateBhpFromThp(const WellState< Scalar, IndexTraits > &well_state, const std::vector< EvalWell > &rates, const Well &well, const SummaryState &summaryState, const Scalar rho, DeferredLogger &deferred_logger) const
std::optional< Scalar > estimateStableBhp(const WellState< Scalar, IndexTraits > &well_state, const Well &well, const std::vector< Scalar > &rates, const Scalar rho, const SummaryState &summaryState) const
int number_of_local_perforations_
Definition: WellInterfaceGeneric.hpp:341
FluidSystem::Scalar wsolvent_
Definition: WellInterfaceGeneric.hpp:381
Definition: WellInterfaceIndices.hpp:34
bool stoppedOrZeroRateTarget(const GroupStateHelperType &groupStateHelper) const
Definition: WellInterface_impl.hpp:1736
bool updateWellOperabilityFromWellEq(const Simulator &simulator, const GroupStateHelperType &groupStateHelper)
Definition: WellInterface_impl.hpp:1288
void checkWellOperability(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper)
Definition: WellInterface_impl.hpp:1160
void updateWellOperability(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper)
Definition: WellInterface_impl.hpp:1250
bool solveWellWithOperabilityCheck(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:601
Scalar woxygen() const
Definition: WellInterface_impl.hpp:166
IndividualOrGroup
Definition: WellInterface.hpp:245
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, WellStateType &well_state, const bool solving_with_zero_rate)
Definition: WellInterface_impl.hpp:974
Scalar computeConnectionDFactor(const int perf, const IntensiveQuantities &intQuants, const SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1983
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
Scalar wfoam() const
Definition: WellInterface_impl.hpp:128
bool updateWellControlAndStatusLocalIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const Scalar WQTotal, WellStateType &well_state, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate)
Definition: WellInterface_impl.hpp:287
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_idx, Callback &extendEval) const
Definition: WellInterface_impl.hpp:2071
std::vector< RateVector > connectionRates_
Definition: WellInterface.hpp:363
bool solveWellForTesting(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:840
void computeConnLevelProdInd(const FluidState &fs, const std::function< Scalar(const Scalar)> &connPICalc, const std::vector< Scalar > &mobility, Scalar *connPI) const
Definition: WellInterface_impl.hpp:2237
void gliftBeginTimeStepWellTestUpdateALQ(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, GLiftEclWells &ecl_well_map, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:1188
Scalar volumetricSurfaceRateForConnection(int cellIdx, int phaseIdx) const
Definition: WellInterface_impl.hpp:1140
virtual void init(const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step)
Definition: WellInterface_impl.hpp:96
std::optional< Scalar > computeBhpAtThpLimitProdWithAlqUsingIPR(const Simulator &simulator, const WellStateType &well_state, Scalar bhp, const SummaryState &summary_state, const Scalar alq_value)
Definition: WellInterface_impl.hpp:2202
void getTw(std::vector< Value > &wi, const int perf, const IntensiveQuantities &intQuants, const Value &trans_mult, const SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1885
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:2084
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
std::vector< Scalar > initialWellRateFractions(const Simulator &ebosSimulator, const WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1747
void solveWellEquation(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:896
void updateConnectionDFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:1963
Eval getPerfCellPressure(const FluidState &fs) const
Definition: WellInterface_impl.hpp:2056
void initializeProducerWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1797
virtual void updateWellStateWithTarget(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1318
void addCellRates(std::map< int, RateVector > &cellRates_) const
Definition: WellInterface_impl.hpp:1121
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
static constexpr bool has_solvent
Definition: WellInterface.hpp:113
bool wellUnderZeroRateTarget(const GroupStateHelperType &groupStateHelper) const
Definition: WellInterface_impl.hpp:1704
GetPropType< TypeTag, Properties::RateVector > RateVector
Definition: WellInterface.hpp:90
void updateConnectionTransmissibilityFactor(const Simulator &simulator, SingleWellStateType &ws) const
Definition: WellInterface_impl.hpp:2025
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:2271
typename GasLiftGroupInfo< Scalar, IndexTraits >::GLiftEclWells GLiftEclWells
Definition: WellInterface.hpp:92
std::unique_ptr< GasLiftSingleWell > initializeGliftWellTest_(const Simulator &simulator, WellStateType &well_state, const GroupState< Scalar > &group_state, GLiftEclWells &ecl_well_map, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:2305
std::optional< Scalar > estimateOperableBhp(const Simulator &ebos_simulator, const double dt, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, WellStateType &well_state)
Definition: WellInterface_impl.hpp:721
Scalar wsalt() const
Definition: WellInterface_impl.hpp:142
bool solveWellWithZeroRate(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:808
bool wellUnderZeroGroupRateTarget(const GroupStateHelperType &groupStateHelper, const std::optional< bool > group_control=std::nullopt) const
Definition: WellInterface_impl.hpp:1721
bool solveWellWithBhp(const Simulator &simulator, const double dt, const Scalar bhp, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:758
void prepareWellBeforeAssembling(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:994
void wellTesting(const Simulator &simulator, const double simulation_time, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, WellTestState &welltest_state, GLiftEclWells &ecl_well_map, std::map< std::string, double > &open_times)
Definition: WellInterface_impl.hpp:406
typename Base::Eval Eval
Definition: WellInterface.hpp:96
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_conservation_quantities, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Constructor.
Definition: WellInterface_impl.hpp:59
bool updateWellStateWithTHPTargetProd(const Simulator &simulator, WellStateType &well_state, const GroupStateHelperType &groupStateHelper) const
Definition: WellInterface_impl.hpp:2170
bool iterateWellEquations(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:532
Scalar wpolymer() const
Definition: WellInterface_impl.hpp:112
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
bool updateWellControl(const Simulator &simulator, const IndividualOrGroup iog, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:190
Scalar wurea() const
Definition: WellInterface_impl.hpp:178
void assembleWellEq(const Simulator &simulator, const double dt, const GroupStateHelperType &groupStateHelper, WellStateType &well_state)
Definition: WellInterface_impl.hpp:958
BlackOilFluidState< Eval, FluidSystem, energyModuleType !=EnergyModules::NoTemperature, energyModuleType==EnergyModules::FullyImplicitThermal, Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, has_solvent, Indices::numPhases > FluidState
Definition: WellInterface.hpp:140
Scalar wmicrobes() const
Definition: WellInterface_impl.hpp:154
virtual void scaleSegmentRatesAndPressure(WellStateType &well_state) const
Definition: WellInterface_impl.hpp:1310
static constexpr bool has_zFraction
Definition: WellInterface.hpp:114
Definition: WellState.hpp:66
constexpr int numPhases() const
The number of phases present.
Definition: WellState.hpp:249
const SingleWellState< Scalar, IndexTraits > & well(std::size_t well_index) const
Definition: WellState.hpp:290
bool isOpen(const std::string &name) const
Definition: WellState.hpp:197
@ NONE
Definition: DeferredLogger.hpp:46
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Phase
Phase indices for reservoir coupling, we currently only support black-oil phases (oil,...
Definition: ReservoirCoupling.hpp:156
Definition: blackoilbioeffectsmodules.hh:45
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Static data associated with a well perforation.
Definition: PerforationData.hpp:30