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