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