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