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