MultisegmentWell_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21// Improve IDE experience
22#ifndef OPM_MULTISEGMENTWELL_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_IMPL_HEADER_INCLUDED
24#include <config.h>
26#endif
27
28#include <opm/common/Exceptions.hpp>
29#include <opm/common/OpmLog/OpmLog.hpp>
30
31#include <opm/input/eclipse/Schedule/MSW/Segment.hpp>
32#include <opm/input/eclipse/Schedule/MSW/Valve.hpp>
33#include <opm/input/eclipse/Schedule/MSW/WellSegments.hpp>
34#include <opm/input/eclipse/Schedule/Well/Connection.hpp>
35#include <opm/input/eclipse/Schedule/Well/WellConnections.hpp>
36
37#include <opm/input/eclipse/Units/Units.hpp>
38
39#include <opm/material/densead/EvaluationFormat.hpp>
40
44
45#include <algorithm>
46#include <cstddef>
47#include <string>
48
49#if COMPILE_BDA_BRIDGE && (HAVE_CUDA || HAVE_OPENCL)
51#endif
52
53namespace Opm
54{
55
56
57 template <typename TypeTag>
59 MultisegmentWell(const Well& well,
60 const ParallelWellInfo<Scalar>& pw_info,
61 const int time_step,
62 const ModelParameters& param,
63 const RateConverterType& rate_converter,
64 const int pvtRegionIdx,
65 const int num_components,
66 const int num_phases,
67 const int index_of_well,
68 const std::vector<PerforationData<Scalar>>& perf_data)
69 : Base(well, pw_info, time_step, param, rate_converter, pvtRegionIdx, num_components, num_phases, index_of_well, perf_data)
70 , MSWEval(static_cast<WellInterfaceIndices<FluidSystem,Indices>&>(*this))
71 , regularize_(false)
72 , segment_fluid_initial_(this->numberOfSegments(), std::vector<Scalar>(this->num_components_, 0.0))
73 {
74 // not handling solvent or polymer for now with multisegment well
75 if constexpr (has_solvent) {
76 OPM_THROW(std::runtime_error, "solvent is not supported by multisegment well yet");
77 }
78
79 if constexpr (has_polymer) {
80 OPM_THROW(std::runtime_error, "polymer is not supported by multisegment well yet");
81 }
82
83 if constexpr (Base::has_energy) {
84 OPM_THROW(std::runtime_error, "energy is not supported by multisegment well yet");
85 }
86
87 if constexpr (Base::has_foam) {
88 OPM_THROW(std::runtime_error, "foam is not supported by multisegment well yet");
89 }
90
91 if constexpr (Base::has_brine) {
92 OPM_THROW(std::runtime_error, "brine is not supported by multisegment well yet");
93 }
94
95 if constexpr (Base::has_watVapor) {
96 OPM_THROW(std::runtime_error, "water evaporation is not supported by multisegment well yet");
97 }
98
99 if(this->rsRvInj() > 0) {
100 OPM_THROW(std::runtime_error,
101 "dissolved gas/ vapporized oil in injected oil/gas not supported by multisegment well yet."
102 " \n See (WCONINJE item 10 / WCONHIST item 8)");
103 }
104 if constexpr (!Indices::oilEnabled && Indices::numPhases > 1) {
105 OPM_THROW(std::runtime_error, "water + gas case not supported by multisegment well yet");
106 }
107
108 this->thp_update_iterations = true;
109 }
110
111
112
113
114
115 template <typename TypeTag>
116 void
118 init(const PhaseUsage* phase_usage_arg,
119 const std::vector<Scalar>& depth_arg,
120 const Scalar gravity_arg,
121 const int num_cells,
122 const std::vector< Scalar >& B_avg,
123 const bool changed_to_open_this_step)
124 {
125 Base::init(phase_usage_arg, depth_arg, gravity_arg, num_cells, B_avg, changed_to_open_this_step);
126
127 // TODO: for StandardWell, we need to update the perf depth here using depth_arg.
128 // for MultisegmentWell, it is much more complicated.
129 // It can be specified directly, it can be calculated from the segment depth,
130 // it can also use the cell center, which is the same for StandardWell.
131 // For the last case, should we update the depth with the depth_arg? For the
132 // future, it can be a source of wrong result with Multisegment well.
133 // An indicator from the opm-parser should indicate what kind of depth we should use here.
134
135 // \Note: we do not update the depth here. And it looks like for now, we only have the option to use
136 // specified perforation depth
137 this->initMatrixAndVectors(num_cells);
138
139 // calculate the depth difference between the perforations and the perforated grid block
140 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
141 const int cell_idx = this->well_cells_[perf];
142 this->cell_perforation_depth_diffs_[perf] = depth_arg[cell_idx] - this->perf_depth_[perf];
143 }
144 }
145
146
147
148
149
150 template <typename TypeTag>
151 void
154 {
155 this->primary_variables_.init();
156 }
157
158
159
160
161
162 template <typename TypeTag>
163 void
165 updatePrimaryVariables(const Simulator& simulator,
166 const WellState<Scalar>& well_state,
167 DeferredLogger& deferred_logger)
168 {
169 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
170 this->primary_variables_.update(well_state, stop_or_zero_rate_target);
171 }
172
173
174
175
176
177
178 template <typename TypeTag>
179 void
182 const GroupState<Scalar>& group_state,
183 WellState<Scalar>& well_state,
184 DeferredLogger& deferred_logger) const
185 {
186 Base::updateWellStateWithTarget(simulator, group_state, well_state, deferred_logger);
187 // scale segment rates based on the wellRates
188 // and segment pressure based on bhp
189 this->scaleSegmentRatesWithWellRates(this->segments_.inlets(),
190 this->segments_.perforations(),
191 well_state);
192 this->scaleSegmentPressuresWithBhp(well_state);
193 }
194
195
196
197
198
199 template <typename TypeTag>
202 getWellConvergence(const Simulator& /* simulator */,
203 const WellState<Scalar>& well_state,
204 const std::vector<Scalar>& B_avg,
205 DeferredLogger& deferred_logger,
206 const bool relax_tolerance) const
207 {
208 return this->MSWEval::getWellConvergence(well_state,
209 B_avg,
210 deferred_logger,
211 this->param_.max_residual_allowed_,
212 this->param_.tolerance_wells_,
213 this->param_.relaxed_tolerance_flow_well_,
214 this->param_.tolerance_pressure_ms_wells_,
215 this->param_.relaxed_tolerance_pressure_ms_well_,
216 relax_tolerance,
217 this->wellIsStopped());
218
219 }
220
221
222
223
224
225 template <typename TypeTag>
226 void
228 apply(const BVector& x, BVector& Ax) const
229 {
230 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
231 return;
232 }
233
234 if (this->param_.matrix_add_well_contributions_) {
235 // Contributions are already in the matrix itself
236 return;
237 }
238
239 this->linSys_.apply(x, Ax);
240 }
241
242
243
244
245
246 template <typename TypeTag>
247 void
249 apply(BVector& r) const
250 {
251 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
252 return;
253 }
254
255 this->linSys_.apply(r);
256 }
257
258
259
260 template <typename TypeTag>
261 void
264 const BVector& x,
265 WellState<Scalar>& well_state,
266 DeferredLogger& deferred_logger)
267 {
268 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) {
269 return;
270 }
271
272 BVectorWell xw(1);
273 this->linSys_.recoverSolutionWell(x, xw);
274 updateWellState(simulator, xw, well_state, deferred_logger);
275 }
276
277
278
279
280
281 template <typename TypeTag>
282 void
284 computeWellPotentials(const Simulator& simulator,
285 const WellState<Scalar>& well_state,
286 std::vector<Scalar>& well_potentials,
287 DeferredLogger& deferred_logger)
288 {
289 const auto [compute_potential, bhp_controlled_well] =
290 this->WellInterfaceGeneric<Scalar>::computeWellPotentials(well_potentials, well_state);
291
292 if (!compute_potential) {
293 return;
294 }
295
296 debug_cost_counter_ = 0;
297 bool converged_implicit = false;
298 if (this->param_.local_well_solver_control_switching_) {
299 converged_implicit = computeWellPotentialsImplicit(simulator, well_potentials, deferred_logger);
300 if (!converged_implicit) {
301 deferred_logger.debug("Implicit potential calculations failed for well "
302 + this->name() + ", reverting to original aproach.");
303 }
304 }
305 if (!converged_implicit) {
306 // does the well have a THP related constraint?
307 const auto& summaryState = simulator.vanguard().summaryState();
308 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
309 computeWellRatesAtBhpLimit(simulator, well_potentials, deferred_logger);
310 } else {
311 well_potentials = computeWellPotentialWithTHP(
312 well_state, simulator, deferred_logger);
313 }
314 }
315 deferred_logger.debug("Cost in iterations of finding well potential for well "
316 + this->name() + ": " + std::to_string(debug_cost_counter_));
317
318 this->checkNegativeWellPotentials(well_potentials,
319 this->param_.check_well_operability_,
320 deferred_logger);
321 }
322
323
324
325
326 template<typename TypeTag>
327 void
330 std::vector<Scalar>& well_flux,
331 DeferredLogger& deferred_logger) const
332 {
333 if (this->well_ecl_.isInjector()) {
334 const auto controls = this->well_ecl_.injectionControls(simulator.vanguard().summaryState());
335 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, well_flux, deferred_logger);
336 } else {
337 const auto controls = this->well_ecl_.productionControls(simulator.vanguard().summaryState());
338 computeWellRatesWithBhpIterations(simulator, controls.bhp_limit, well_flux, deferred_logger);
339 }
340 }
341
342 template<typename TypeTag>
343 void
345 computeWellRatesWithBhp(const Simulator& simulator,
346 const Scalar& bhp,
347 std::vector<Scalar>& well_flux,
348 DeferredLogger& deferred_logger) const
349 {
350 const int np = this->number_of_phases_;
351
352 well_flux.resize(np, 0.0);
353 const bool allow_cf = this->getAllowCrossFlow();
354 const int nseg = this->numberOfSegments();
355 const WellState<Scalar>& well_state = simulator.problem().wellModel().wellState();
356 const auto& ws = well_state.well(this->indexOfWell());
357 auto segments_copy = ws.segments;
358 segments_copy.scale_pressure(bhp);
359 const auto& segment_pressure = segments_copy.pressure;
360 for (int seg = 0; seg < nseg; ++seg) {
361 for (const int perf : this->segments_.perforations()[seg]) {
362 const int cell_idx = this->well_cells_[perf];
363 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
364 // flux for each perforation
365 std::vector<Scalar> mob(this->num_components_, 0.);
366 getMobility(simulator, perf, mob, deferred_logger);
367 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
368 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
369 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
370 const Scalar seg_pressure = segment_pressure[seg];
371 std::vector<Scalar> cq_s(this->num_components_, 0.);
372 Scalar perf_press = 0.0;
373 PerforationRates<Scalar> perf_rates;
374 computePerfRate(intQuants, mob, Tw, seg, perf, seg_pressure,
375 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
376
377 for(int p = 0; p < np; ++p) {
378 well_flux[this->modelCompIdxToFlowCompIdx(p)] += cq_s[p];
379 }
380 }
381 }
382 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
383 }
384
385
386 template<typename TypeTag>
387 void
390 const Scalar& bhp,
391 std::vector<Scalar>& well_flux,
392 DeferredLogger& deferred_logger) const
393 {
394 // creating a copy of the well itself, to avoid messing up the explicit information
395 // during this copy, the only information not copied properly is the well controls
396 MultisegmentWell<TypeTag> well_copy(*this);
397 well_copy.debug_cost_counter_ = 0;
398
399 // store a copy of the well state, we don't want to update the real well state
400 WellState<Scalar> well_state_copy = simulator.problem().wellModel().wellState();
401 const auto& group_state = simulator.problem().wellModel().groupState();
402 auto& ws = well_state_copy.well(this->index_of_well_);
403
404 // Get the current controls.
405 const auto& summary_state = simulator.vanguard().summaryState();
406 auto inj_controls = well_copy.well_ecl_.isInjector()
407 ? well_copy.well_ecl_.injectionControls(summary_state)
408 : Well::InjectionControls(0);
409 auto prod_controls = well_copy.well_ecl_.isProducer()
410 ? well_copy.well_ecl_.productionControls(summary_state) :
411 Well::ProductionControls(0);
412
413 // Set current control to bhp, and bhp value in state, modify bhp limit in control object.
414 if (well_copy.well_ecl_.isInjector()) {
415 inj_controls.bhp_limit = bhp;
416 ws.injection_cmode = Well::InjectorCMode::BHP;
417 } else {
418 prod_controls.bhp_limit = bhp;
419 ws.production_cmode = Well::ProducerCMode::BHP;
420 }
421 ws.bhp = bhp;
422 well_copy.scaleSegmentPressuresWithBhp(well_state_copy);
423
424 // initialized the well rates with the potentials i.e. the well rates based on bhp
425 const int np = this->number_of_phases_;
426 bool trivial = true;
427 for (int phase = 0; phase < np; ++phase){
428 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
429 }
430 if (!trivial) {
431 const Scalar sign = well_copy.well_ecl_.isInjector() ? 1.0 : -1.0;
432 for (int phase = 0; phase < np; ++phase) {
433 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
434 }
435 }
436 well_copy.scaleSegmentRatesWithWellRates(this->segments_.inlets(),
437 this->segments_.perforations(),
438 well_state_copy);
439
440 well_copy.calculateExplicitQuantities(simulator, well_state_copy, deferred_logger);
441 const double dt = simulator.timeStepSize();
442 // iterate to get a solution at the given bhp.
443 well_copy.iterateWellEqWithControl(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state,
444 deferred_logger);
445
446 // compute the potential and store in the flux vector.
447 well_flux.clear();
448 well_flux.resize(np, 0.0);
449 for (int compIdx = 0; compIdx < this->num_components_; ++compIdx) {
450 const EvalWell rate = well_copy.primary_variables_.getQs(compIdx);
451 well_flux[this->modelCompIdxToFlowCompIdx(compIdx)] = rate.value();
452 }
453 debug_cost_counter_ += well_copy.debug_cost_counter_;
454 }
455
456
457
458 template<typename TypeTag>
459 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
462 const Simulator& simulator,
463 DeferredLogger& deferred_logger) const
464 {
465 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
466 const auto& summary_state = simulator.vanguard().summaryState();
467
468 const auto& well = this->well_ecl_;
469 if (well.isInjector()) {
470 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, summary_state, deferred_logger);
471 if (bhp_at_thp_limit) {
472 const auto& controls = well.injectionControls(summary_state);
473 const Scalar bhp = std::min(*bhp_at_thp_limit,
474 static_cast<Scalar>(controls.bhp_limit));
475 computeWellRatesWithBhpIterations(simulator, bhp, potentials, deferred_logger);
476 deferred_logger.debug("Converged thp based potential calculation for well "
477 + this->name() + ", at bhp = " + std::to_string(bhp));
478 } else {
479 deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL",
480 "Failed in getting converged thp based potential calculation for well "
481 + this->name() + ". Instead the bhp based value is used");
482 const auto& controls = well.injectionControls(summary_state);
483 const Scalar bhp = controls.bhp_limit;
484 computeWellRatesWithBhpIterations(simulator, bhp, potentials, deferred_logger);
485 }
486 } else {
487 auto bhp_at_thp_limit = computeBhpAtThpLimitProd(
488 well_state, simulator, summary_state, deferred_logger);
489 if (bhp_at_thp_limit) {
490 const auto& controls = well.productionControls(summary_state);
491 const Scalar bhp = std::max(*bhp_at_thp_limit,
492 static_cast<Scalar>(controls.bhp_limit));
493 computeWellRatesWithBhpIterations(simulator, bhp, potentials, deferred_logger);
494 deferred_logger.debug("Converged thp based potential calculation for well "
495 + this->name() + ", at bhp = " + std::to_string(bhp));
496 } else {
497 deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL",
498 "Failed in getting converged thp based potential calculation for well "
499 + this->name() + ". Instead the bhp based value is used");
500 const auto& controls = well.productionControls(summary_state);
501 const Scalar bhp = controls.bhp_limit;
502 computeWellRatesWithBhpIterations(simulator, bhp, potentials, deferred_logger);
503 }
504 }
505
506 return potentials;
507 }
508
509 template<typename TypeTag>
510 bool
513 std::vector<Scalar>& well_potentials,
514 DeferredLogger& deferred_logger) const
515 {
516 // Create a copy of the well.
517 // TODO: check if we can avoid taking multiple copies. Call from updateWellPotentials
518 // is allready a copy, but not from other calls.
519 MultisegmentWell<TypeTag> well_copy(*this);
520 well_copy.debug_cost_counter_ = 0;
521
522 // store a copy of the well state, we don't want to update the real well state
523 WellState<Scalar> well_state_copy = simulator.problem().wellModel().wellState();
524 const auto& group_state = simulator.problem().wellModel().groupState();
525 auto& ws = well_state_copy.well(this->index_of_well_);
526
527 // get current controls
528 const auto& summary_state = simulator.vanguard().summaryState();
529 auto inj_controls = well_copy.well_ecl_.isInjector()
530 ? well_copy.well_ecl_.injectionControls(summary_state)
531 : Well::InjectionControls(0);
532 auto prod_controls = well_copy.well_ecl_.isProducer()
533 ? well_copy.well_ecl_.productionControls(summary_state)
534 : Well::ProductionControls(0);
535
536 // prepare/modify well state and control
537 well_copy.prepareForPotentialCalculations(summary_state, well_state_copy, inj_controls, prod_controls);
538
539 well_copy.scaleSegmentPressuresWithBhp(well_state_copy);
540
541 // initialize rates from previous potentials
542 const int np = this->number_of_phases_;
543 bool trivial = true;
544 for (int phase = 0; phase < np; ++phase){
545 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
546 }
547 if (!trivial) {
548 const Scalar sign = well_copy.well_ecl_.isInjector() ? 1.0 : -1.0;
549 for (int phase = 0; phase < np; ++phase) {
550 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
551 }
552 }
553 well_copy.scaleSegmentRatesWithWellRates(this->segments_.inlets(),
554 this->segments_.perforations(),
555 well_state_copy);
556
557 well_copy.calculateExplicitQuantities(simulator, well_state_copy, deferred_logger);
558 const double dt = simulator.timeStepSize();
559 // solve equations
560 bool converged = false;
561 if (this->well_ecl_.isProducer() && this->wellHasTHPConstraints(summary_state)) {
562 converged = well_copy.solveWellWithTHPConstraint(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
563 } else {
564 converged = well_copy.iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
565 }
566
567 // fetch potentials (sign is updated on the outside).
568 well_potentials.clear();
569 well_potentials.resize(np, 0.0);
570 for (int compIdx = 0; compIdx < this->num_components_; ++compIdx) {
571 const EvalWell rate = well_copy.primary_variables_.getQs(compIdx);
572 well_potentials[this->modelCompIdxToFlowCompIdx(compIdx)] = rate.value();
573 }
574 debug_cost_counter_ += well_copy.debug_cost_counter_;
575 return converged;
576 }
577
578 template <typename TypeTag>
579 void
582 WellState<Scalar>& well_state,
583 DeferredLogger& deferred_logger)
584 {
585 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
586
587 // We assemble the well equations, then we check the convergence,
588 // which is why we do not put the assembleWellEq here.
589 try{
590 const BVectorWell dx_well = this->linSys_.solve();
591
592 updateWellState(simulator, dx_well, well_state, deferred_logger);
593 }
594 catch(const NumericalProblem& exp) {
595 // Add information about the well and log to deferred logger
596 // (Logging done inside of solve() method will only be seen if
597 // this is the process with rank zero)
598 deferred_logger.problem("In MultisegmentWell::solveEqAndUpdateWellState for well "
599 + this->name() +": "+exp.what());
600 throw;
601 }
602 }
603
604
605
606
607
608 template <typename TypeTag>
609 void
612 {
613 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
614
615 std::vector<Scalar> kr(this->number_of_phases_, 0.0);
616 std::vector<Scalar> density(this->number_of_phases_, 0.0);
617
618 const int cell_idx = this->well_cells_[perf];
619 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
620 const auto& fs = intQuants.fluidState();
621
622 Scalar sum_kr = 0.;
623
624 const PhaseUsage& pu = this->phaseUsage();
625 if (pu.phase_used[Water]) {
626 const int water_pos = pu.phase_pos[Water];
627 kr[water_pos] = intQuants.relativePermeability(FluidSystem::waterPhaseIdx).value();
628 sum_kr += kr[water_pos];
629 density[water_pos] = fs.density(FluidSystem::waterPhaseIdx).value();
630 }
631
632 if (pu.phase_used[Oil]) {
633 const int oil_pos = pu.phase_pos[Oil];
634 kr[oil_pos] = intQuants.relativePermeability(FluidSystem::oilPhaseIdx).value();
635 sum_kr += kr[oil_pos];
636 density[oil_pos] = fs.density(FluidSystem::oilPhaseIdx).value();
637 }
638
639 if (pu.phase_used[Gas]) {
640 const int gas_pos = pu.phase_pos[Gas];
641 kr[gas_pos] = intQuants.relativePermeability(FluidSystem::gasPhaseIdx).value();
642 sum_kr += kr[gas_pos];
643 density[gas_pos] = fs.density(FluidSystem::gasPhaseIdx).value();
644 }
645
646 assert(sum_kr != 0.);
647
648 // calculate the average density
649 Scalar average_density = 0.;
650 for (int p = 0; p < this->number_of_phases_; ++p) {
651 average_density += kr[p] * density[p];
652 }
653 average_density /= sum_kr;
654
655 this->cell_perforation_pressure_diffs_[perf] = this->gravity_ * average_density * this->cell_perforation_depth_diffs_[perf];
656 }
657 }
658
659
660
661
662
663 template <typename TypeTag>
664 void
667 {
668 for (int seg = 0; seg < this->numberOfSegments(); ++seg) {
669 // TODO: trying to reduce the times for the surfaceVolumeFraction calculation
670 const Scalar surface_volume = getSegmentSurfaceVolume(simulator, seg).value();
671 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
672 segment_fluid_initial_[seg][comp_idx] = surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx).value();
673 }
674 }
675 }
676
677
678
679
680
681 template <typename TypeTag>
682 void
684 updateWellState(const Simulator& simulator,
685 const BVectorWell& dwells,
686 WellState<Scalar>& well_state,
687 DeferredLogger& deferred_logger,
688 const Scalar relaxation_factor)
689 {
690 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
691
692 const Scalar dFLimit = this->param_.dwell_fraction_max_;
693 const Scalar max_pressure_change = this->param_.max_pressure_change_ms_wells_;
694 const bool stop_or_zero_rate_target =
695 this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
696 this->primary_variables_.updateNewton(dwells,
697 relaxation_factor,
698 dFLimit,
699 stop_or_zero_rate_target,
700 max_pressure_change);
701
702 const auto& summary_state = simulator.vanguard().summaryState();
703 this->primary_variables_.copyToWellState(*this, getRefDensity(),
704 stop_or_zero_rate_target,
705 well_state,
706 summary_state,
707 deferred_logger);
708
709 {
710 auto& ws = well_state.well(this->index_of_well_);
711 this->segments_.copyPhaseDensities(ws.pu, ws.segments);
712 }
713
714 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.well(this->index_of_well_));
715 }
716
717
718
719
720
721 template <typename TypeTag>
722 void
725 const WellState<Scalar>& well_state,
726 DeferredLogger& deferred_logger)
727 {
728 updatePrimaryVariables(simulator, well_state, deferred_logger);
729 initPrimaryVariablesEvaluation();
730 computePerfCellPressDiffs(simulator);
731 computeInitialSegmentFluids(simulator);
732 }
733
734
735
736
737
738 template<typename TypeTag>
739 void
741 updateProductivityIndex(const Simulator& simulator,
742 const WellProdIndexCalculator<Scalar>& wellPICalc,
743 WellState<Scalar>& well_state,
744 DeferredLogger& deferred_logger) const
745 {
746 auto fluidState = [&simulator, this](const int perf)
747 {
748 const auto cell_idx = this->well_cells_[perf];
749 return simulator.model()
750 .intensiveQuantities(cell_idx, /*timeIdx=*/ 0).fluidState();
751 };
752
753 const int np = this->number_of_phases_;
754 auto setToZero = [np](Scalar* x) -> void
755 {
756 std::fill_n(x, np, 0.0);
757 };
758
759 auto addVector = [np](const Scalar* src, Scalar* dest) -> void
760 {
761 std::transform(src, src + np, dest, dest, std::plus<>{});
762 };
763
764 auto& ws = well_state.well(this->index_of_well_);
765 auto& perf_data = ws.perf_data;
766 auto* connPI = perf_data.prod_index.data();
767 auto* wellPI = ws.productivity_index.data();
768
769 setToZero(wellPI);
770
771 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
772 auto subsetPerfID = 0;
773
774 for ( const auto& perf : *this->perf_data_){
775 auto allPerfID = perf.ecl_index;
776
777 auto connPICalc = [&wellPICalc, allPerfID](const Scalar mobility) -> Scalar
778 {
779 return wellPICalc.connectionProdIndStandard(allPerfID, mobility);
780 };
781
782 std::vector<Scalar> mob(this->num_components_, 0.0);
783 getMobility(simulator, static_cast<int>(subsetPerfID), mob, deferred_logger);
784
785 const auto& fs = fluidState(subsetPerfID);
786 setToZero(connPI);
787
788 if (this->isInjector()) {
789 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
790 mob, connPI, deferred_logger);
791 }
792 else { // Production or zero flow rate
793 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
794 }
795
796 addVector(connPI, wellPI);
797
798 ++subsetPerfID;
799 connPI += np;
800 }
801
802 assert (static_cast<int>(subsetPerfID) == this->number_of_perforations_ &&
803 "Internal logic error in processing connections for PI/II");
804 }
805
806
807
808
809
810 template<typename TypeTag>
813 connectionDensity(const int globalConnIdx,
814 [[maybe_unused]] const int openConnIdx) const
815 {
816 // Simple approximation: Mixture density at reservoir connection is
817 // mixture density at connection's segment.
818
819 const auto segNum = this->wellEcl()
820 .getConnections()[globalConnIdx].segment();
821
822 const auto segIdx = this->wellEcl()
823 .getSegments().segmentNumberToIndex(segNum);
824
825 return this->segments_.density(segIdx).value();
826 }
827
828
829
830
831
832 template<typename TypeTag>
833 void
836 {
837 this->linSys_.extract(jacobian);
838 }
839
840
841 template<typename TypeTag>
842 void
845 const BVector& weights,
846 const int pressureVarIndex,
847 const bool use_well_weights,
848 const WellState<Scalar>& well_state) const
849 {
850 // Add the pressure contribution to the cpr system for the well
851 this->linSys_.extractCPRPressureMatrix(jacobian,
852 weights,
853 pressureVarIndex,
854 use_well_weights,
855 *this,
856 this->SPres,
857 well_state);
858 }
859
860
861 template<typename TypeTag>
862 template<class Value>
863 void
865 computePerfRate(const Value& pressure_cell,
866 const Value& rs,
867 const Value& rv,
868 const std::vector<Value>& b_perfcells,
869 const std::vector<Value>& mob_perfcells,
870 const std::vector<Scalar>& Tw,
871 const int perf,
872 const Value& segment_pressure,
873 const Value& segment_density,
874 const bool& allow_cf,
875 const std::vector<Value>& cmix_s,
876 std::vector<Value>& cq_s,
877 Value& perf_press,
878 PerforationRates<Scalar>& perf_rates,
879 DeferredLogger& deferred_logger) const
880 {
881 // pressure difference between the segment and the perforation
882 const Value perf_seg_press_diff = this->gravity() * segment_density *
883 this->segments_.perforation_depth_diff(perf);
884 // pressure difference between the perforation and the grid cell
885 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[perf];
886
887 // perforation pressure is the wellbore pressure corrected to perforation depth
888 // (positive sign due to convention in segments_.perforation_depth_diff() )
889 perf_press = segment_pressure + perf_seg_press_diff;
890
891 // cell pressure corrected to perforation depth
892 const Value cell_press_at_perf = pressure_cell - cell_perf_press_diff;
893
894 // Pressure drawdown (also used to determine direction of flow)
895 const Value drawdown = cell_press_at_perf - perf_press;
896
897 // producing perforations
898 if (drawdown > 0.0) {
899 // Do nothing if crossflow is not allowed
900 if (!allow_cf && this->isInjector()) {
901 return;
902 }
903
904 // compute component volumetric rates at standard conditions
905 for (int comp_idx = 0; comp_idx < this->numComponents(); ++comp_idx) {
906 const Value cq_p = - Tw[comp_idx] * (mob_perfcells[comp_idx] * drawdown);
907 cq_s[comp_idx] = b_perfcells[comp_idx] * cq_p;
908 }
909
910 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
911 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
912 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
913 const Value cq_s_oil = cq_s[oilCompIdx];
914 const Value cq_s_gas = cq_s[gasCompIdx];
915 cq_s[gasCompIdx] += rs * cq_s_oil;
916 cq_s[oilCompIdx] += rv * cq_s_gas;
917 }
918 } else { // injecting perforations
919 // Do nothing if crossflow is not allowed
920 if (!allow_cf && this->isProducer()) {
921 return;
922 }
923
924 // for injecting perforations, we use total mobility
925 Value total_mob = mob_perfcells[0];
926 for (int comp_idx = 1; comp_idx < this->numComponents(); ++comp_idx) {
927 total_mob += mob_perfcells[comp_idx];
928 }
929
930 // compute volume ratio between connection and at standard conditions
931 Value volume_ratio = 0.0;
932 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
933 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
934 volume_ratio += cmix_s[waterCompIdx] / b_perfcells[waterCompIdx];
935 }
936
937 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
938 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
939 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
940
941 // Incorporate RS/RV factors if both oil and gas active
942 // TODO: not sure we use rs rv from the perforation cells when handling injecting perforations
943 // basically, for injecting perforations, the wellbore is the upstreaming side.
944 const Value d = 1.0 - rv * rs;
945
946 if (getValue(d) == 0.0) {
947 OPM_DEFLOG_PROBLEM(NumericalProblem,
948 fmt::format("Zero d value obtained for well {} "
949 "during flux calculation with rs {} and rv {}",
950 this->name(), rs, rv),
951 deferred_logger);
952 }
953
954 const Value tmp_oil = (cmix_s[oilCompIdx] - rv * cmix_s[gasCompIdx]) / d;
955 volume_ratio += tmp_oil / b_perfcells[oilCompIdx];
956
957 const Value tmp_gas = (cmix_s[gasCompIdx] - rs * cmix_s[oilCompIdx]) / d;
958 volume_ratio += tmp_gas / b_perfcells[gasCompIdx];
959 } else { // not having gas and oil at the same time
960 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
961 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
962 volume_ratio += cmix_s[oilCompIdx] / b_perfcells[oilCompIdx];
963 }
964 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
965 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
966 volume_ratio += cmix_s[gasCompIdx] / b_perfcells[gasCompIdx];
967 }
968 }
969 // injecting connections total volumerates at standard conditions
970 for (int componentIdx = 0; componentIdx < this->numComponents(); ++componentIdx) {
971 const Value cqt_i = - Tw[componentIdx] * (total_mob * drawdown);
972 Value cqt_is = cqt_i / volume_ratio;
973 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
974 }
975 } // end for injection perforations
976
977 // calculating the perforation solution gas rate and solution oil rates
978 if (this->isProducer()) {
979 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
980 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
981 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
982 // TODO: the formulations here remain to be tested with cases with strong crossflow through production wells
983 // s means standard condition, r means reservoir condition
984 // q_os = q_or * b_o + rv * q_gr * b_g
985 // q_gs = q_gr * g_g + rs * q_or * b_o
986 // d = 1.0 - rs * rv
987 // q_or = 1 / (b_o * d) * (q_os - rv * q_gs)
988 // q_gr = 1 / (b_g * d) * (q_gs - rs * q_os)
989
990 const Scalar d = 1.0 - getValue(rv) * getValue(rs);
991 // vaporized oil into gas
992 // rv * q_gr * b_g = rv * (q_gs - rs * q_os) / d
993 perf_rates.vap_oil = getValue(rv) * (getValue(cq_s[gasCompIdx]) - getValue(rs) * getValue(cq_s[oilCompIdx])) / d;
994 // dissolved of gas in oil
995 // rs * q_or * b_o = rs * (q_os - rv * q_gs) / d
996 perf_rates.dis_gas = getValue(rs) * (getValue(cq_s[oilCompIdx]) - getValue(rv) * getValue(cq_s[gasCompIdx])) / d;
997 }
998 }
999 }
1000
1001 template <typename TypeTag>
1002 template<class Value>
1003 void
1005 computePerfRate(const IntensiveQuantities& int_quants,
1006 const std::vector<Value>& mob_perfcells,
1007 const std::vector<Scalar>& Tw,
1008 const int seg,
1009 const int perf,
1010 const Value& segment_pressure,
1011 const bool& allow_cf,
1012 std::vector<Value>& cq_s,
1013 Value& perf_press,
1014 PerforationRates<Scalar>& perf_rates,
1015 DeferredLogger& deferred_logger) const
1016
1017 {
1018 auto obtain = [this](const Eval& value)
1019 {
1020 if constexpr (std::is_same_v<Value, Scalar>) {
1021 static_cast<void>(this); // suppress clang warning
1022 return getValue(value);
1023 } else {
1024 return this->extendEval(value);
1025 }
1026 };
1027 auto obtainN = [](const auto& value)
1028 {
1029 if constexpr (std::is_same_v<Value, Scalar>) {
1030 return getValue(value);
1031 } else {
1032 return value;
1033 }
1034 };
1035 const auto& fs = int_quants.fluidState();
1036
1037 const Value pressure_cell = obtain(this->getPerfCellPressure(fs));
1038 const Value rs = obtain(fs.Rs());
1039 const Value rv = obtain(fs.Rv());
1040
1041 // not using number_of_phases_ because of solvent
1042 std::vector<Value> b_perfcells(this->num_components_, 0.0);
1043
1044 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
1045 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1046 continue;
1047 }
1048
1049 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
1050 b_perfcells[compIdx] = obtain(fs.invB(phaseIdx));
1051 }
1052
1053 std::vector<Value> cmix_s(this->numComponents(), 0.0);
1054 for (int comp_idx = 0; comp_idx < this->numComponents(); ++comp_idx) {
1055 cmix_s[comp_idx] = obtainN(this->primary_variables_.surfaceVolumeFraction(seg, comp_idx));
1056 }
1057
1058 this->computePerfRate(pressure_cell,
1059 rs,
1060 rv,
1061 b_perfcells,
1062 mob_perfcells,
1063 Tw,
1064 perf,
1065 segment_pressure,
1066 obtainN(this->segments_.density(seg)),
1067 allow_cf,
1068 cmix_s,
1069 cq_s,
1070 perf_press,
1071 perf_rates,
1072 deferred_logger);
1073 }
1074
1075 template <typename TypeTag>
1076 void
1078 computeSegmentFluidProperties(const Simulator& simulator, DeferredLogger& deferred_logger)
1079 {
1080 // TODO: the concept of phases and components are rather confusing in this function.
1081 // needs to be addressed sooner or later.
1082
1083 // get the temperature for later use. It is only useful when we are not handling
1084 // thermal related simulation
1085 // basically, it is a single value for all the segments
1086
1087 EvalWell temperature;
1088 EvalWell saltConcentration;
1089 // not sure how to handle the pvt region related to segment
1090 // for the current approach, we use the pvt region of the first perforated cell
1091 // although there are some text indicating using the pvt region of the lowest
1092 // perforated cell
1093 // TODO: later to investigate how to handle the pvt region
1094 int pvt_region_index;
1095 {
1096 // using the first perforated cell
1097 const int cell_idx = this->well_cells_[0];
1098 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
1099 const auto& fs = intQuants.fluidState();
1100 temperature.setValue(fs.temperature(FluidSystem::oilPhaseIdx).value());
1101 saltConcentration = this->extendEval(fs.saltConcentration());
1102 pvt_region_index = fs.pvtRegionIndex();
1103 }
1104
1105 this->segments_.computeFluidProperties(temperature,
1106 saltConcentration,
1107 this->primary_variables_,
1108 pvt_region_index,
1109 deferred_logger);
1110 }
1111
1112 template <typename TypeTag>
1113 template<class Value>
1114 void
1116 getMobility(const Simulator& simulator,
1117 const int perf,
1118 std::vector<Value>& mob,
1119 DeferredLogger& deferred_logger) const
1120 {
1121 auto obtain = [this](const Eval& value)
1122 {
1123 if constexpr (std::is_same_v<Value, Scalar>) {
1124 static_cast<void>(this); // suppress clang warning
1125 return getValue(value);
1126 } else {
1127 return this->extendEval(value);
1128 }
1129 };
1130
1131 WellInterface<TypeTag>::getMobility(simulator, perf, mob, obtain, deferred_logger);
1132
1133 if (this->isInjector() && this->well_ecl_.getInjMultMode() != Well::InjMultMode::NONE) {
1134 const auto perf_ecl_index = this->perforationData()[perf].ecl_index;
1135 const Connection& con = this->well_ecl_.getConnections()[perf_ecl_index];
1136 const int seg = this->segmentNumberToIndex(con.segment());
1137 // from the reference results, it looks like MSW uses segment pressure instead of BHP here
1138 // Note: this is against the documented definition.
1139 // we can change this depending on what we want
1140 const Scalar segment_pres = this->primary_variables_.getSegmentPressure(seg).value();
1141 const Scalar perf_seg_press_diff = this->gravity() * this->segments_.density(seg).value()
1142 * this->segments_.perforation_depth_diff(perf);
1143 const Scalar perf_press = segment_pres + perf_seg_press_diff;
1144 const Scalar multiplier = this->getInjMult(perf, segment_pres, perf_press);
1145 for (std::size_t i = 0; i < mob.size(); ++i) {
1146 mob[i] *= multiplier;
1147 }
1148 }
1149 }
1150
1151
1152
1153 template<typename TypeTag>
1156 getRefDensity() const
1157 {
1158 return this->segments_.getRefDensity();
1159 }
1160
1161 template<typename TypeTag>
1162 void
1165 const Simulator& simulator,
1166 DeferredLogger& deferred_logger)
1167 {
1168 const auto& summaryState = simulator.vanguard().summaryState();
1169 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1170 // Crude but works: default is one atmosphere.
1171 // TODO: a better way to detect whether the BHP is defaulted or not
1172 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
1173 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
1174 // if the BHP limit is not defaulted or the well does not have a THP limit
1175 // we need to check the BHP limit
1176 Scalar total_ipr_mass_rate = 0.0;
1177 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1178 {
1179 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1180 continue;
1181 }
1182
1183 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
1184 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
1185
1186 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
1187 total_ipr_mass_rate += ipr_rate * rho;
1188 }
1189 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
1190 this->operability_status_.operable_under_only_bhp_limit = false;
1191 }
1192
1193 // checking whether running under BHP limit will violate THP limit
1194 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
1195 // option 1: calculate well rates based on the BHP limit.
1196 // option 2: stick with the above IPR curve
1197 // we use IPR here
1198 std::vector<Scalar> well_rates_bhp_limit;
1199 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1200
1201 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1202 const Scalar thp = WellBhpThpCalculator(*this).calculateThpFromBhp(well_rates_bhp_limit,
1203 bhp_limit,
1204 this->getRefDensity(),
1205 this->wellEcl().alq_value(summaryState),
1206 thp_limit,
1207 deferred_logger);
1208 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1209 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
1210 }
1211 }
1212 } else {
1213 // defaulted BHP and there is a THP constraint
1214 // default BHP limit is about 1 atm.
1215 // when applied the hydrostatic pressure correction dp,
1216 // most likely we get a negative value (bhp + dp)to search in the VFP table,
1217 // which is not desirable.
1218 // we assume we can operate under defaulted BHP limit and will violate the THP limit
1219 // when operating under defaulted BHP limit.
1220 this->operability_status_.operable_under_only_bhp_limit = true;
1221 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
1222 }
1223 }
1224
1225
1226
1227 template<typename TypeTag>
1228 void
1230 updateIPR(const Simulator& simulator, DeferredLogger& deferred_logger) const
1231 {
1232 // TODO: not handling solvent related here for now
1233
1234 // initialize all the values to be zero to begin with
1235 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
1236 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
1237
1238 const int nseg = this->numberOfSegments();
1239 std::vector<Scalar> seg_dp(nseg, 0.0);
1240 for (int seg = 0; seg < nseg; ++seg) {
1241 // calculating the perforation rate for each perforation that belongs to this segment
1242 const Scalar dp = this->getSegmentDp(seg,
1243 this->segments_.density(seg).value(),
1244 seg_dp);
1245 seg_dp[seg] = dp;
1246 for (const int perf : this->segments_.perforations()[seg]) {
1247 std::vector<Scalar> mob(this->num_components_, 0.0);
1248
1249 // TODO: maybe we should store the mobility somewhere, so that we only need to calculate it one per iteration
1250 getMobility(simulator, perf, mob, deferred_logger);
1251
1252 const int cell_idx = this->well_cells_[perf];
1253 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1254 const auto& fs = int_quantities.fluidState();
1255 // pressure difference between the segment and the perforation
1256 const Scalar perf_seg_press_diff = this->segments_.getPressureDiffSegPerf(seg, perf);
1257 // pressure difference between the perforation and the grid cell
1258 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[perf];
1259 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
1260
1261 // calculating the b for the connection
1262 std::vector<Scalar> b_perf(this->num_components_);
1263 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
1264 if (!FluidSystem::phaseIsActive(phase)) {
1265 continue;
1266 }
1267 const unsigned comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phase));
1268 b_perf[comp_idx] = fs.invB(phase).value();
1269 }
1270
1271 // the pressure difference between the connection and BHP
1272 const Scalar h_perf = cell_perf_press_diff + perf_seg_press_diff + dp;
1273 const Scalar pressure_diff = pressure_cell - h_perf;
1274
1275 // do not take into consideration the crossflow here.
1276 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
1277 deferred_logger.debug("CROSSFLOW_IPR",
1278 "cross flow found when updateIPR for well " + this->name());
1279 }
1280
1281 // the well index associated with the connection
1282 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quantities, cell_idx);
1283 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1284 const std::vector<Scalar> tw_perf = this->wellIndex(perf, int_quantities, trans_mult, wellstate_nupcol);
1285 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
1286 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
1287 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1288 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
1289 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
1290 ipr_b_perf[comp_idx] += tw_mob;
1291 }
1292
1293 // we need to handle the rs and rv when both oil and gas are present
1294 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1295 const unsigned oil_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
1296 const unsigned gas_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
1297 const Scalar rs = (fs.Rs()).value();
1298 const Scalar rv = (fs.Rv()).value();
1299
1300 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
1301 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
1302
1303 ipr_a_perf[gas_comp_idx] += dis_gas_a;
1304 ipr_a_perf[oil_comp_idx] += vap_oil_a;
1305
1306 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
1307 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
1308
1309 ipr_b_perf[gas_comp_idx] += dis_gas_b;
1310 ipr_b_perf[oil_comp_idx] += vap_oil_b;
1311 }
1312
1313 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
1314 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
1315 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
1316 }
1317 }
1318 }
1319 }
1320
1321 template<typename TypeTag>
1322 void
1324 updateIPRImplicit(const Simulator& simulator,
1325 WellState<Scalar>& well_state,
1326 DeferredLogger& deferred_logger)
1327 {
1328 // Compute IPR based on *converged* well-equation:
1329 // For a component rate r the derivative dr/dbhp is obtained by
1330 // dr/dbhp = - (partial r/partial x) * inv(partial Eq/partial x) * (partial Eq/partial bhp_target)
1331 // where Eq(x)=0 is the well equation setup with bhp control and primary variables x
1332
1333 // We shouldn't have zero rates at this stage, but check
1334 bool zero_rates;
1335 auto rates = well_state.well(this->index_of_well_).surface_rates;
1336 zero_rates = true;
1337 for (std::size_t p = 0; p < rates.size(); ++p) {
1338 zero_rates &= rates[p] == 0.0;
1339 }
1340 auto& ws = well_state.well(this->index_of_well_);
1341 if (zero_rates) {
1342 const auto msg = fmt::format("updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
1343 deferred_logger.debug(msg);
1344 /*
1345 // could revert to standard approach here:
1346 updateIPR(simulator, deferred_logger);
1347 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx){
1348 const int idx = this->modelCompIdxToFlowCompIdx(comp_idx);
1349 ws.implicit_ipr_a[idx] = this->ipr_a_[comp_idx];
1350 ws.implicit_ipr_b[idx] = this->ipr_b_[comp_idx];
1351 }
1352 return;
1353 */
1354 }
1355 const auto& group_state = simulator.problem().wellModel().groupState();
1356
1357 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
1358 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
1359 //WellState well_state_copy = well_state;
1360 auto inj_controls = Well::InjectionControls(0);
1361 auto prod_controls = Well::ProductionControls(0);
1362 prod_controls.addControl(Well::ProducerCMode::BHP);
1363 prod_controls.bhp_limit = well_state.well(this->index_of_well_).bhp;
1364
1365 // Set current control to bhp, and bhp value in state, modify bhp limit in control object.
1366 const auto cmode = ws.production_cmode;
1367 ws.production_cmode = Well::ProducerCMode::BHP;
1368 const double dt = simulator.timeStepSize();
1369 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
1370
1371 BVectorWell rhs(this->numberOfSegments());
1372 rhs = 0.0;
1373 rhs[0][SPres] = -1.0;
1374
1375 const BVectorWell x_well = this->linSys_.solve(rhs);
1376 constexpr int num_eq = MSWEval::numWellEq;
1377 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx){
1378 const EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
1379 const int idx = this->modelCompIdxToFlowCompIdx(comp_idx);
1380 for (size_t pvIdx = 0; pvIdx < num_eq; ++pvIdx) {
1381 // well primary variable derivatives in EvalWell start at position Indices::numEq
1382 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
1383 }
1384 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
1385 }
1386 // reset cmode
1387 ws.production_cmode = cmode;
1388 }
1389
1390 template<typename TypeTag>
1391 void
1394 const WellState<Scalar>& well_state,
1395 DeferredLogger& deferred_logger)
1396 {
1397 const auto& summaryState = simulator.vanguard().summaryState();
1398 const auto obtain_bhp = this->isProducer()
1399 ? computeBhpAtThpLimitProd(
1400 well_state, simulator, summaryState, deferred_logger)
1401 : computeBhpAtThpLimitInj(simulator, summaryState, deferred_logger);
1402
1403 if (obtain_bhp) {
1404 this->operability_status_.can_obtain_bhp_with_thp_limit = true;
1405
1406 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1407 this->operability_status_.obey_bhp_limit_with_thp_limit = (*obtain_bhp >= bhp_limit);
1408
1409 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1410 if (this->isProducer() && *obtain_bhp < thp_limit) {
1411 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1412 + " bars is SMALLER than thp limit "
1413 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1414 + " bars as a producer for well " + this->name();
1415 deferred_logger.debug(msg);
1416 }
1417 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1418 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1419 + " bars is LARGER than thp limit "
1420 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1421 + " bars as a injector for well " + this->name();
1422 deferred_logger.debug(msg);
1423 }
1424 } else {
1425 // Shutting wells that can not find bhp value from thp
1426 // when under THP control
1427 this->operability_status_.can_obtain_bhp_with_thp_limit = false;
1428 this->operability_status_.obey_bhp_limit_with_thp_limit = false;
1429 if (!this->wellIsStopped()) {
1430 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1431 deferred_logger.debug(" could not find bhp value at thp limit "
1432 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1433 + " bar for well " + this->name() + ", the well might need to be closed ");
1434 }
1435 }
1436 }
1437
1438
1439
1440
1441
1442 template<typename TypeTag>
1443 bool
1445 iterateWellEqWithControl(const Simulator& simulator,
1446 const double dt,
1447 const Well::InjectionControls& inj_controls,
1448 const Well::ProductionControls& prod_controls,
1449 WellState<Scalar>& well_state,
1450 const GroupState<Scalar>& group_state,
1451 DeferredLogger& deferred_logger)
1452 {
1453 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return true;
1454
1455 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1456
1457 {
1458 // getWellFiniteResiduals returns false for nan/inf residuals
1459 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1460 if(!isFinite)
1461 return false;
1462 }
1463
1464 std::vector<std::vector<Scalar> > residual_history;
1465 std::vector<Scalar> measure_history;
1466 int it = 0;
1467 // relaxation factor
1468 Scalar relaxation_factor = 1.;
1469 const Scalar min_relaxation_factor = 0.6;
1470 bool converged = false;
1471 int stagnate_count = 0;
1472 bool relax_convergence = false;
1473 this->regularize_ = false;
1474 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1475
1476 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
1477
1478 BVectorWell dx_well;
1479 try{
1480 dx_well = this->linSys_.solve();
1481 }
1482 catch(const NumericalProblem& exp) {
1483 // Add information about the well and log to deferred logger
1484 // (Logging done inside of solve() method will only be seen if
1485 // this is the process with rank zero)
1486 deferred_logger.problem("In MultisegmentWell::iterateWellEqWithControl for well "
1487 + this->name() +": "+exp.what());
1488 throw;
1489 }
1490
1491 if (it > this->param_.strict_inner_iter_wells_) {
1492 relax_convergence = true;
1493 this->regularize_ = true;
1494 }
1495
1496 const auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
1497 if (report.converged()) {
1498 converged = true;
1499 break;
1500 }
1501
1502 {
1503 // getFinteWellResiduals returns false for nan/inf residuals
1504 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1505 if (!isFinite)
1506 return false;
1507
1508 residual_history.push_back(residuals);
1509 measure_history.push_back(this->getResidualMeasureValue(well_state,
1510 residual_history[it],
1511 this->param_.tolerance_wells_,
1512 this->param_.tolerance_pressure_ms_wells_,
1513 deferred_logger) );
1514 }
1515
1516
1517 bool is_oscillate = false;
1518 bool is_stagnate = false;
1519
1520 this->detectOscillations(measure_history, is_oscillate, is_stagnate);
1521 // TODO: maybe we should have more sophisticated strategy to recover the relaxation factor,
1522 // for example, to recover it to be bigger
1523
1524 if (is_oscillate || is_stagnate) {
1525 // HACK!
1526 std::ostringstream sstr;
1527 if (relaxation_factor == min_relaxation_factor) {
1528 // Still stagnating, terminate iterations if 5 iterations pass.
1529 ++stagnate_count;
1530 if (stagnate_count == 6) {
1531 sstr << " well " << this->name() << " observes severe stagnation and/or oscillation. We relax the tolerance and check for convergence. \n";
1532 const auto reportStag = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, true);
1533 if (reportStag.converged()) {
1534 converged = true;
1535 sstr << " well " << this->name() << " manages to get converged with relaxed tolerances in " << it << " inner iterations";
1536 deferred_logger.debug(sstr.str());
1537 return converged;
1538 }
1539 }
1540 }
1541
1542 // a factor value to reduce the relaxation_factor
1543 const Scalar reduction_mutliplier = 0.9;
1544 relaxation_factor = std::max(relaxation_factor * reduction_mutliplier, min_relaxation_factor);
1545
1546 // debug output
1547 if (is_stagnate) {
1548 sstr << " well " << this->name() << " observes stagnation in inner iteration " << it << "\n";
1549
1550 }
1551 if (is_oscillate) {
1552 sstr << " well " << this->name() << " observes oscillation in inner iteration " << it << "\n";
1553 }
1554 sstr << " relaxation_factor is " << relaxation_factor << " now\n";
1555
1556 this->regularize_ = true;
1557 deferred_logger.debug(sstr.str());
1558 }
1559 updateWellState(simulator, dx_well, well_state, deferred_logger, relaxation_factor);
1560 initPrimaryVariablesEvaluation();
1561 }
1562
1563 // TODO: we should decide whether to keep the updated well_state, or recover to use the old well_state
1564 if (converged) {
1565 std::ostringstream sstr;
1566 sstr << " Well " << this->name() << " converged in " << it << " inner iterations.";
1567 if (relax_convergence)
1568 sstr << " (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ << " iterations)";
1569 deferred_logger.debug(sstr.str());
1570 } else {
1571 std::ostringstream sstr;
1572 sstr << " Well " << this->name() << " did not converge in " << it << " inner iterations.";
1573#define EXTRA_DEBUG_MSW 0
1574#if EXTRA_DEBUG_MSW
1575 sstr << "***** Outputting the residual history for well " << this->name() << " during inner iterations:";
1576 for (int i = 0; i < it; ++i) {
1577 const auto& residual = residual_history[i];
1578 sstr << " residual at " << i << "th iteration ";
1579 for (const auto& res : residual) {
1580 sstr << " " << res;
1581 }
1582 sstr << " " << measure_history[i] << " \n";
1583 }
1584#endif
1585#undef EXTRA_DEBUG_MSW
1586 deferred_logger.debug(sstr.str());
1587 }
1588
1589 return converged;
1590 }
1591
1592
1593 template<typename TypeTag>
1594 bool
1596 iterateWellEqWithSwitching(const Simulator& simulator,
1597 const double dt,
1598 const Well::InjectionControls& inj_controls,
1599 const Well::ProductionControls& prod_controls,
1600 WellState<Scalar>& well_state,
1601 const GroupState<Scalar>& group_state,
1602 DeferredLogger& deferred_logger,
1603 const bool fixed_control /*false*/,
1604 const bool fixed_status /*false*/)
1605 {
1606 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1607
1608 {
1609 // getWellFiniteResiduals returns false for nan/inf residuals
1610 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1611 if(!isFinite)
1612 return false;
1613 }
1614
1615 std::vector<std::vector<Scalar> > residual_history;
1616 std::vector<Scalar> measure_history;
1617 int it = 0;
1618 // relaxation factor
1619 Scalar relaxation_factor = 1.;
1620 const Scalar min_relaxation_factor = 0.6;
1621 bool converged = false;
1622 [[maybe_unused]] int stagnate_count = 0;
1623 bool relax_convergence = false;
1624 this->regularize_ = false;
1625 const auto& summary_state = simulator.vanguard().summaryState();
1626
1627 // Always take a few (more than one) iterations after a switch before allowing a new switch
1628 // The optimal number here is subject to further investigation, but it has been observerved
1629 // that unless this number is >1, we may get stuck in a cycle
1630 const int min_its_after_switch = 3;
1631 int its_since_last_switch = min_its_after_switch;
1632 int switch_count= 0;
1633 // if we fail to solve eqs, we reset status/operability before leaving
1634 const auto well_status_orig = this->wellStatus_;
1635 const auto operability_orig = this->operability_status_;
1636 auto well_status_cur = well_status_orig;
1637 // don't allow opening wells that are stopped from schedule or has a stopped well state
1638 const bool allow_open = this->well_ecl_.getStatus() == WellStatus::OPEN &&
1639 well_state.well(this->index_of_well_).status == WellStatus::OPEN;
1640 // don't allow switcing for wells under zero rate target or requested fixed status and control
1641 const bool allow_switching = !this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) &&
1642 (!fixed_control || !fixed_status) && allow_open;
1643 bool changed = false;
1644 bool final_check = false;
1645 // well needs to be set operable or else solving/updating of re-opened wells is skipped
1646 this->operability_status_.resetOperability();
1647 this->operability_status_.solvable = true;
1648
1649 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1650 its_since_last_switch++;
1651 if (allow_switching && its_since_last_switch >= min_its_after_switch){
1652 const Scalar wqTotal = this->primary_variables_.getWQTotal().value();
1653 changed = this->updateWellControlAndStatusLocalIteration(simulator, well_state, group_state,
1654 inj_controls, prod_controls, wqTotal,
1655 deferred_logger, fixed_control, fixed_status);
1656 if (changed){
1657 its_since_last_switch = 0;
1658 switch_count++;
1659 if (well_status_cur != this->wellStatus_) {
1660 well_status_cur = this->wellStatus_;
1661 }
1662 }
1663 if (!changed && final_check) {
1664 break;
1665 } else {
1666 final_check = false;
1667 }
1668 }
1669
1670 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
1671
1672 const BVectorWell dx_well = this->linSys_.solve();
1673
1674 if (it > this->param_.strict_inner_iter_wells_) {
1675 relax_convergence = true;
1676 this->regularize_ = true;
1677 }
1678
1679 const auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
1680 converged = report.converged();
1681 if (converged) {
1682 // if equations are sufficiently linear they might converge in less than min_its_after_switch
1683 // in this case, make sure all constraints are satisfied before returning
1684 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
1685 final_check = true;
1686 its_since_last_switch = min_its_after_switch;
1687 } else {
1688 break;
1689 }
1690 }
1691
1692 // getFinteWellResiduals returns false for nan/inf residuals
1693 {
1694 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1695 if (!isFinite)
1696 return false;
1697
1698 residual_history.push_back(residuals);
1699 }
1700
1701 if (!converged) {
1702 measure_history.push_back(this->getResidualMeasureValue(well_state,
1703 residual_history[it],
1704 this->param_.tolerance_wells_,
1705 this->param_.tolerance_pressure_ms_wells_,
1706 deferred_logger));
1707
1708 bool is_oscillate = false;
1709 bool is_stagnate = false;
1710
1711 this->detectOscillations(measure_history, is_oscillate, is_stagnate);
1712 // TODO: maybe we should have more sophisticated strategy to recover the relaxation factor,
1713 // for example, to recover it to be bigger
1714
1715 if (is_oscillate || is_stagnate) {
1716 // HACK!
1717 std::string message;
1718 if (relaxation_factor == min_relaxation_factor) {
1719 ++stagnate_count;
1720 if (false) { // this disables the usage of the relaxed tolerance
1721 fmt::format_to(std::back_inserter(message), " Well {} observes severe stagnation and/or oscillation."
1722 " We relax the tolerance and check for convergence. \n", this->name());
1723 const auto reportStag = getWellConvergence(simulator, well_state, Base::B_avg_,
1724 deferred_logger, true);
1725 if (reportStag.converged()) {
1726 converged = true;
1727 fmt::format_to(std::back_inserter(message), " Well {} manages to get converged with relaxed tolerances in {} inner iterations", this->name(), it);
1728 deferred_logger.debug(message);
1729 return converged;
1730 }
1731 }
1732 }
1733
1734 // a factor value to reduce the relaxation_factor
1735 constexpr Scalar reduction_mutliplier = 0.9;
1736 relaxation_factor = std::max(relaxation_factor * reduction_mutliplier, min_relaxation_factor);
1737
1738 // debug output
1739 if (is_stagnate) {
1740 fmt::format_to(std::back_inserter(message), " well {} observes stagnation in inner iteration {}\n", this->name(), it);
1741 }
1742 if (is_oscillate) {
1743 fmt::format_to(std::back_inserter(message), " well {} observes oscillation in inner iteration {}\n", this->name(), it);
1744 }
1745 fmt::format_to(std::back_inserter(message), " relaxation_factor is {} now\n", relaxation_factor);
1746
1747 this->regularize_ = true;
1748 deferred_logger.debug(message);
1749 }
1750 }
1751 updateWellState(simulator, dx_well, well_state, deferred_logger, relaxation_factor);
1752 initPrimaryVariablesEvaluation();
1753 }
1754
1755 if (converged) {
1756 if (allow_switching){
1757 // update operability if status change
1758 const bool is_stopped = this->wellIsStopped();
1759 if (this->wellHasTHPConstraints(summary_state)){
1760 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
1761 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
1762 } else {
1763 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
1764 }
1765 }
1766 std::string message = fmt::format(" Well {} converged in {} inner iterations ("
1767 "{} control/status switches).", this->name(), it, switch_count);
1768 if (relax_convergence) {
1769 message.append(fmt::format(" (A relaxed tolerance was used after {} iterations)",
1770 this->param_.strict_inner_iter_wells_));
1771 }
1772 deferred_logger.debug(message);
1773 } else {
1774 this->wellStatus_ = well_status_orig;
1775 this->operability_status_ = operability_orig;
1776 const std::string message = fmt::format(" Well {} did not converge in {} inner iterations ("
1777 "{} control/status switches).", this->name(), it, switch_count);
1778 deferred_logger.debug(message);
1779 this->primary_variables_.outputLowLimitPressureSegments(deferred_logger);
1780 }
1781
1782 return converged;
1783 }
1784
1785
1786 template<typename TypeTag>
1787 void
1790 const double dt,
1791 const Well::InjectionControls& inj_controls,
1792 const Well::ProductionControls& prod_controls,
1793 WellState<Scalar>& well_state,
1794 const GroupState<Scalar>& group_state,
1795 DeferredLogger& deferred_logger)
1796 {
1797 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1798
1799 // update the upwinding segments
1800 this->segments_.updateUpwindingSegments(this->primary_variables_);
1801
1802 // calculate the fluid properties needed.
1803 computeSegmentFluidProperties(simulator, deferred_logger);
1804
1805 // clear all entries
1806 this->linSys_.clear();
1807
1808 auto& ws = well_state.well(this->index_of_well_);
1809 ws.phase_mixing_rates.fill(0.0);
1810
1811 // for the black oil cases, there will be four equations,
1812 // the first three of them are the mass balance equations, the last one is the pressure equations.
1813 //
1814 // but for the top segment, the pressure equation will be the well control equation, and the other three will be the same.
1815
1816 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1817
1818 const int nseg = this->numberOfSegments();
1819
1820 for (int seg = 0; seg < nseg; ++seg) {
1821 // calculating the accumulation term
1822 // TODO: without considering the efficiency factor for now
1823 {
1824 const EvalWell segment_surface_volume = getSegmentSurfaceVolume(simulator, seg);
1825
1826 // Add a regularization_factor to increase the accumulation term
1827 // This will make the system less stiff and help convergence for
1828 // difficult cases
1829 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
1830 // for each component
1831 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1832 const EvalWell accumulation_term = regularization_factor * (segment_surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx)
1833 - segment_fluid_initial_[seg][comp_idx]) / dt;
1835 assembleAccumulationTerm(seg, comp_idx, accumulation_term, this->linSys_);
1836 }
1837 }
1838 // considering the contributions due to flowing out from the segment
1839 {
1840 const int seg_upwind = this->segments_.upwinding_segment(seg);
1841 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1842 const EvalWell segment_rate =
1843 this->primary_variables_.getSegmentRateUpwinding(seg,
1844 seg_upwind,
1845 comp_idx) *
1846 this->well_efficiency_factor_;
1848 assembleOutflowTerm(seg, seg_upwind, comp_idx, segment_rate, this->linSys_);
1849 }
1850 }
1851
1852 // considering the contributions from the inlet segments
1853 {
1854 for (const int inlet : this->segments_.inlets()[seg]) {
1855 const int inlet_upwind = this->segments_.upwinding_segment(inlet);
1856 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1857 const EvalWell inlet_rate =
1858 this->primary_variables_.getSegmentRateUpwinding(inlet,
1859 inlet_upwind,
1860 comp_idx) *
1861 this->well_efficiency_factor_;
1863 assembleInflowTerm(seg, inlet, inlet_upwind, comp_idx, inlet_rate, this->linSys_);
1864 }
1865 }
1866 }
1867
1868 // calculating the perforation rate for each perforation that belongs to this segment
1869 const EvalWell seg_pressure = this->primary_variables_.getSegmentPressure(seg);
1870 auto& perf_data = ws.perf_data;
1871 auto& perf_rates = perf_data.phase_rates;
1872 auto& perf_press_state = perf_data.pressure;
1873 for (const int perf : this->segments_.perforations()[seg]) {
1874 const int cell_idx = this->well_cells_[perf];
1875 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1876 std::vector<EvalWell> mob(this->num_components_, 0.0);
1877 getMobility(simulator, perf, mob, deferred_logger);
1878 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quants, cell_idx);
1879 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1880 const std::vector<Scalar> Tw = this->wellIndex(perf, int_quants, trans_mult, wellstate_nupcol);
1881 std::vector<EvalWell> cq_s(this->num_components_, 0.0);
1882 EvalWell perf_press;
1883 PerforationRates<Scalar> perfRates;
1884 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
1885 allow_cf, cq_s, perf_press, perfRates, deferred_logger);
1886
1887 // updating the solution gas rate and solution oil rate
1888 if (this->isProducer()) {
1889 ws.phase_mixing_rates[ws.dissolved_gas] += perfRates.dis_gas;
1890 ws.phase_mixing_rates[ws.vaporized_oil] += perfRates.vap_oil;
1891 perf_data.phase_mixing_rates[perf][ws.dissolved_gas] = perfRates.dis_gas;
1892 perf_data.phase_mixing_rates[perf][ws.vaporized_oil] = perfRates.vap_oil;
1893 }
1894
1895 // store the perf pressure and rates
1896 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1897 perf_rates[perf*this->number_of_phases_ + this->modelCompIdxToFlowCompIdx(comp_idx)] = cq_s[comp_idx].value();
1898 }
1899 perf_press_state[perf] = perf_press.value();
1900
1901 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1902 // the cq_s entering mass balance equations need to consider the efficiency factors.
1903 const EvalWell cq_s_effective = cq_s[comp_idx] * this->well_efficiency_factor_;
1904
1905 this->connectionRates_[perf][comp_idx] = Base::restrictEval(cq_s_effective);
1906
1908 assemblePerforationEq(seg, cell_idx, comp_idx, cq_s_effective, this->linSys_);
1909 }
1910 }
1911
1912 // the fourth equation, the pressure drop equation
1913 if (seg == 0) { // top segment, pressure equation is the control equation
1914 const auto& summaryState = simulator.vanguard().summaryState();
1915 const Schedule& schedule = simulator.vanguard().schedule();
1916 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1918 assembleControlEq(well_state,
1919 group_state,
1920 schedule,
1921 summaryState,
1922 inj_controls,
1923 prod_controls,
1924 getRefDensity(),
1925 this->primary_variables_,
1926 this->linSys_,
1927 stopped_or_zero_target,
1928 deferred_logger);
1929 } else {
1930 const UnitSystem& unit_system = simulator.vanguard().eclState().getDeckUnitSystem();
1931 const auto& summary_state = simulator.vanguard().summaryState();
1932 this->assemblePressureEq(seg, unit_system, well_state, summary_state, this->param_.use_average_density_ms_wells_, deferred_logger);
1933 }
1934 }
1935
1936 this->linSys_.createSolver();
1937 }
1938
1939
1940
1941
1942 template<typename TypeTag>
1943 bool
1945 openCrossFlowAvoidSingularity(const Simulator& simulator) const
1946 {
1947 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1948 }
1949
1950
1951 template<typename TypeTag>
1952 bool
1954 allDrawDownWrongDirection(const Simulator& simulator) const
1955 {
1956 bool all_drawdown_wrong_direction = true;
1957 const int nseg = this->numberOfSegments();
1958
1959 for (int seg = 0; seg < nseg; ++seg) {
1960 const EvalWell segment_pressure = this->primary_variables_.getSegmentPressure(seg);
1961 for (const int perf : this->segments_.perforations()[seg]) {
1962
1963 const int cell_idx = this->well_cells_[perf];
1964 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1965 const auto& fs = intQuants.fluidState();
1966
1967 // pressure difference between the segment and the perforation
1968 const EvalWell perf_seg_press_diff = this->segments_.getPressureDiffSegPerf(seg, perf);
1969 // pressure difference between the perforation and the grid cell
1970 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[perf];
1971
1972 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
1973 const Scalar perf_press = pressure_cell - cell_perf_press_diff;
1974 // Pressure drawdown (also used to determine direction of flow)
1975 // TODO: not 100% sure about the sign of the seg_perf_press_diff
1976 const EvalWell drawdown = perf_press - (segment_pressure + perf_seg_press_diff);
1977
1978 // for now, if there is one perforation can produce/inject in the correct
1979 // direction, we consider this well can still produce/inject.
1980 // TODO: it can be more complicated than this to cause wrong-signed rates
1981 if ( (drawdown < 0. && this->isInjector()) ||
1982 (drawdown > 0. && this->isProducer()) ) {
1983 all_drawdown_wrong_direction = false;
1984 break;
1985 }
1986 }
1987 }
1988
1989 return all_drawdown_wrong_direction;
1990 }
1991
1992
1993
1994
1995 template<typename TypeTag>
1996 void
1998 updateWaterThroughput(const double /*dt*/, WellState<Scalar>& /*well_state*/) const
1999 {
2000 }
2001
2002
2003
2004
2005
2006 template<typename TypeTag>
2009 getSegmentSurfaceVolume(const Simulator& simulator, const int seg_idx) const
2010 {
2011 EvalWell temperature;
2012 EvalWell saltConcentration;
2013 int pvt_region_index;
2014 {
2015 // using the pvt region of first perforated cell
2016 // TODO: it should be a member of the WellInterface, initialized properly
2017 const int cell_idx = this->well_cells_[0];
2018 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
2019 const auto& fs = intQuants.fluidState();
2020 temperature.setValue(fs.temperature(FluidSystem::oilPhaseIdx).value());
2021 saltConcentration = this->extendEval(fs.saltConcentration());
2022 pvt_region_index = fs.pvtRegionIndex();
2023 }
2024
2025 return this->segments_.getSurfaceVolume(temperature,
2026 saltConcentration,
2027 this->primary_variables_,
2028 pvt_region_index,
2029 seg_idx);
2030 }
2031
2032
2033 template<typename TypeTag>
2034 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2037 const Simulator& simulator,
2038 const SummaryState& summary_state,
2039 DeferredLogger& deferred_logger) const
2040 {
2042 simulator,
2043 summary_state,
2044 this->getALQ(well_state),
2045 deferred_logger);
2046 }
2047
2048
2049
2050 template<typename TypeTag>
2051 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2054 const SummaryState& summary_state,
2055 const Scalar alq_value,
2056 DeferredLogger& deferred_logger) const
2057 {
2058 // Make the frates() function.
2059 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2060 // Not solving the well equations here, which means we are
2061 // calculating at the current Fg/Fw values of the
2062 // well. This does not matter unless the well is
2063 // crossflowing, and then it is likely still a good
2064 // approximation.
2065 std::vector<Scalar> rates(3);
2066 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2067 return rates;
2068 };
2069
2070 auto bhpAtLimit = WellBhpThpCalculator(*this).
2071 computeBhpAtThpLimitProd(frates,
2072 summary_state,
2073 this->maxPerfPress(simulator),
2074 this->getRefDensity(),
2075 alq_value,
2076 this->getTHPConstraint(summary_state),
2077 deferred_logger);
2078
2079 if (bhpAtLimit)
2080 return bhpAtLimit;
2081
2082 auto fratesIter = [this, &simulator, &deferred_logger](const Scalar bhp) {
2083 // Solver the well iterations to see if we are
2084 // able to get a solution with an update
2085 // solution
2086 std::vector<Scalar> rates(3);
2087 computeWellRatesWithBhpIterations(simulator, bhp, rates, deferred_logger);
2088 return rates;
2089 };
2090
2091 return WellBhpThpCalculator(*this).
2092 computeBhpAtThpLimitProd(fratesIter,
2093 summary_state,
2094 this->maxPerfPress(simulator),
2095 this->getRefDensity(),
2096 alq_value,
2097 this->getTHPConstraint(summary_state),
2098 deferred_logger);
2099 }
2100
2101 template<typename TypeTag>
2102 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2104 computeBhpAtThpLimitInj(const Simulator& simulator,
2105 const SummaryState& summary_state,
2106 DeferredLogger& deferred_logger) const
2107 {
2108 // Make the frates() function.
2109 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2110 // Not solving the well equations here, which means we are
2111 // calculating at the current Fg/Fw values of the
2112 // well. This does not matter unless the well is
2113 // crossflowing, and then it is likely still a good
2114 // approximation.
2115 std::vector<Scalar> rates(3);
2116 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2117 return rates;
2118 };
2119
2120 auto bhpAtLimit = WellBhpThpCalculator(*this).
2121 computeBhpAtThpLimitInj(frates,
2122 summary_state,
2123 this->getRefDensity(),
2124 0.05,
2125 100,
2126 false,
2127 deferred_logger);
2128
2129 if (bhpAtLimit)
2130 return bhpAtLimit;
2131
2132 auto fratesIter = [this, &simulator, &deferred_logger](const Scalar bhp) {
2133 // Solver the well iterations to see if we are
2134 // able to get a solution with an update
2135 // solution
2136 std::vector<Scalar> rates(3);
2137 computeWellRatesWithBhpIterations(simulator, bhp, rates, deferred_logger);
2138 return rates;
2139 };
2140
2141 return WellBhpThpCalculator(*this).
2142 computeBhpAtThpLimitInj(fratesIter,
2143 summary_state,
2144 this->getRefDensity(),
2145 0.05,
2146 100,
2147 false,
2148 deferred_logger);
2149 }
2150
2151
2152
2153
2154
2155 template<typename TypeTag>
2158 maxPerfPress(const Simulator& simulator) const
2159 {
2160 Scalar max_pressure = 0.0;
2161 const int nseg = this->numberOfSegments();
2162 for (int seg = 0; seg < nseg; ++seg) {
2163 for (const int perf : this->segments_.perforations()[seg]) {
2164 const int cell_idx = this->well_cells_[perf];
2165 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2166 const auto& fs = int_quants.fluidState();
2167 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2168 max_pressure = std::max(max_pressure, pressure_cell);
2169 }
2170 }
2171 return max_pressure;
2172 }
2173
2174
2175
2176
2177
2178 template<typename TypeTag>
2179 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2181 computeCurrentWellRates(const Simulator& simulator,
2182 DeferredLogger& deferred_logger) const
2183 {
2184 // Calculate the rates that follow from the current primary variables.
2185 std::vector<Scalar> well_q_s(this->num_components_, 0.0);
2186 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2187 const int nseg = this->numberOfSegments();
2188 for (int seg = 0; seg < nseg; ++seg) {
2189 // calculating the perforation rate for each perforation that belongs to this segment
2190 const Scalar seg_pressure = getValue(this->primary_variables_.getSegmentPressure(seg));
2191 for (const int perf : this->segments_.perforations()[seg]) {
2192 const int cell_idx = this->well_cells_[perf];
2193 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2194 std::vector<Scalar> mob(this->num_components_, 0.0);
2195 getMobility(simulator, perf, mob, deferred_logger);
2196 const Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quants, cell_idx);
2197 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2198 const std::vector<Scalar> Tw = this->wellIndex(perf, int_quants, trans_mult, wellstate_nupcol);
2199 std::vector<Scalar> cq_s(this->num_components_, 0.0);
2200 Scalar perf_press = 0.0;
2201 PerforationRates<Scalar> perf_rates;
2202 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
2203 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
2204 for (int comp = 0; comp < this->num_components_; ++comp) {
2205 well_q_s[comp] += cq_s[comp];
2206 }
2207 }
2208 }
2209 return well_q_s;
2210 }
2211
2212
2213 template <typename TypeTag>
2214 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2216 getPrimaryVars() const
2217 {
2218 const int num_seg = this->numberOfSegments();
2219 constexpr int num_eq = MSWEval::numWellEq;
2220 std::vector<Scalar> retval(num_seg * num_eq);
2221 for (int ii = 0; ii < num_seg; ++ii) {
2222 const auto& pv = this->primary_variables_.value(ii);
2223 std::copy(pv.begin(), pv.end(), retval.begin() + ii * num_eq);
2224 }
2225 return retval;
2226 }
2227
2228
2229
2230
2231 template <typename TypeTag>
2232 int
2234 setPrimaryVars(typename std::vector<Scalar>::const_iterator it)
2235 {
2236 const int num_seg = this->numberOfSegments();
2237 constexpr int num_eq = MSWEval::numWellEq;
2238 std::array<Scalar, num_eq> tmp;
2239 for (int ii = 0; ii < num_seg; ++ii) {
2240 const auto start = it + num_seg * num_eq;
2241 std::copy(start, start + num_eq, tmp.begin());
2242 this->primary_variables_.setValue(ii, tmp);
2243 }
2244 return num_seg * num_eq;
2245 }
2246
2247} // namespace Opm
#define OPM_DEFLOG_PROBLEM(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:61
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
void warning(const std::string &tag, const std::string &message)
void problem(const std::string &tag, const std::string &message)
void debug(const std::string &tag, const std::string &message)
Definition: GroupState.hpp:38
Class handling assemble of the equation system for MultisegmentWell.
Definition: MultisegmentWellAssemble.hpp:44
PrimaryVariables primary_variables_
The primary variables.
Definition: MultisegmentWellEval.hpp:139
void scaleSegmentPressuresWithBhp(WellState< Scalar > &well_state) const
void scaleSegmentRatesWithWellRates(const std::vector< std::vector< int > > &segment_inlets, const std::vector< std::vector< int > > &segment_perforations, WellState< Scalar > &well_state) const
Definition: MultisegmentWell.hpp:36
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2053
std::vector< Scalar > computeWellPotentialWithTHP(const WellState< Scalar > &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:461
void computeWellPotentials(const Simulator &simulator, const WellState< Scalar > &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: MultisegmentWell_impl.hpp:284
void calculateExplicitQuantities(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:724
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:813
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:835
std::optional< Scalar > computeBhpAtThpLimitProd(const WellState< Scalar > &well_state, const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2036
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1156
bool iterateWellEqWithSwitching(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, const bool fixed_control=false, const bool fixed_status=false) override
Definition: MultisegmentWell_impl.hpp:1596
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2216
void solveEqAndUpdateWellState(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:581
ConvergenceReport getWellConvergence(const Simulator &simulator, const WellState< Scalar > &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition: MultisegmentWell_impl.hpp:202
void updatePrimaryVariables(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:165
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:741
void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:389
std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2181
void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: MultisegmentWell_impl.hpp:228
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1116
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellState< Scalar > &well_state, DeferredLogger &deferred_logger, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:684
void computePerfRate(const IntensiveQuantities &int_quants, const std::vector< Value > &mob_perfcells, const std::vector< Scalar > &Tw, const int seg, const int perf, const Value &segment_pressure, const bool &allow_cf, std::vector< Value > &cq_s, Value &perf_press, PerforationRates< Scalar > &perf_rates, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1005
bool iterateWellEqWithControl(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) override
Definition: MultisegmentWell_impl.hpp:1445
void checkOperabilityUnderBHPLimit(const WellState< Scalar > &well_state, const Simulator &ebos_simulator, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1164
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:1945
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1078
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:263
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: MultisegmentWell_impl.hpp:2234
void computeWellRatesWithBhp(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:345
EvalWell getSegmentSurfaceVolume(const Simulator &simulator, const int seg_idx) const
Definition: MultisegmentWell_impl.hpp:2009
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:1954
void updateWellStateWithTarget(const Simulator &simulator, const GroupState< Scalar > &group_state, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
updating the well state based the current control mode
Definition: MultisegmentWell_impl.hpp:181
int debug_cost_counter_
Definition: MultisegmentWell.hpp:174
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellState< Scalar > &well_state) const override
Definition: MultisegmentWell_impl.hpp:844
static constexpr bool has_polymer
Definition: WellInterface.hpp:110
void assembleWellEqWithoutIteration(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) override
Definition: MultisegmentWell_impl.hpp:1789
static constexpr bool has_solvent
Definition: WellInterface.hpp:108
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1393
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:611
void updateIPR(const Simulator &ebos_simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:1230
void updateWaterThroughput(const double dt, WellState< Scalar > &well_state) const override
Definition: MultisegmentWell_impl.hpp:1998
void initPrimaryVariablesEvaluation() override
Definition: MultisegmentWell_impl.hpp:153
bool computeWellPotentialsImplicit(const Simulator &simulator, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:512
void computeWellRatesAtBhpLimit(const Simulator &simulator, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:329
MultisegmentWell(const Well &well, const ParallelWellInfo< Scalar > &pw_info, const int time_step, const ModelParameters &param, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_components, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Definition: MultisegmentWell_impl.hpp:59
Scalar maxPerfPress(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2158
void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1324
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) override
Definition: MultisegmentWell_impl.hpp:118
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2104
void computeInitialSegmentFluids(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:666
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
Class for computing BHP limits.
Definition: WellBhpThpCalculator.hpp:41
Scalar calculateThpFromBhp(const std::vector< Scalar > &rates, const Scalar bhp, const Scalar rho, const std::optional< Scalar > &alq, const Scalar thp_limit, DeferredLogger &deferred_logger) const
Calculates THP from BHP.
Scalar mostStrictBhpFromBhpLimits(const SummaryState &summaryState) const
Obtain the most strict BHP from BHP limits.
void prepareForPotentialCalculations(const SummaryState &summary_state, WellState< Scalar > &well_state, Well::InjectionControls &inj_controls, Well::ProductionControls &prod_controls) const
std::pair< bool, bool > computeWellPotentials(std::vector< Scalar > &well_potentials, const WellState< Scalar > &well_state)
Well well_ecl_
Definition: WellInterfaceGeneric.hpp:273
Definition: WellInterfaceIndices.hpp:34
Definition: WellInterface.hpp:73
static constexpr bool has_brine
Definition: WellInterface.hpp:116
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:78
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1773
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:100
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:97
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:82
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
static constexpr bool has_watVapor
Definition: WellInterface.hpp:117
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:80
static constexpr bool has_foam
Definition: WellInterface.hpp:115
typename Base::Eval Eval
Definition: WellInterface.hpp:95
static constexpr bool has_energy
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:81
bool solveWellWithTHPConstraint(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellState< Scalar > &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger)
Definition: WellInterface_impl.hpp:514
bool thp_update_iterations
Definition: WellInterface.hpp:392
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:84
Definition: WellProdIndexCalculator.hpp:37
Scalar connectionProdIndStandard(const std::size_t connIdx, const Scalar connMobility) const
Definition: WellState.hpp:62
const SingleWellState< Scalar > & well(std::size_t well_index) const
Definition: WellState.hpp:307
@ NONE
Definition: DeferredLogger.hpp:46
void detectOscillations(const std::vector< std::vector< Scalar > > &residualHistory, const int it, const int numPhases, const Scalar relaxRelTol, const int minimumOscillatingPhases, bool &oscillate, bool &stagnate)
Detect oscillation or stagnation in a given residual history.
Definition: blackoilboundaryratevector.hh:37
PhaseUsage phaseUsage(const Phases &phases)
Determine the active phases.
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:147
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: PerforationData.hpp:41
Scalar dis_gas
Definition: PerforationData.hpp:42
Scalar vap_oil
Definition: PerforationData.hpp:44
Definition: BlackoilPhases.hpp:46
std::array< int, MaxNumPhases+NumCryptoPhases > phase_used
Definition: BlackoilPhases.hpp:51
std::array< int, MaxNumPhases+NumCryptoPhases > phase_pos
Definition: BlackoilPhases.hpp:52