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 );
598 }
599
600 // fetch potentials (sign is updated on the outside).
601 well_potentials.clear();
602 well_potentials.resize(np, 0.0);
603 for (int compIdx = 0; compIdx < this->num_conservation_quantities_; ++compIdx) {
604 const EvalWell rate = well_copy.primary_variables_.getQs(compIdx);
605 well_potentials[FluidSystem::activeCompToActivePhaseIdx(compIdx)] = rate.value();
606 }
607 debug_cost_counter_ += well_copy.debug_cost_counter_;
608 return converged;
609 }
610
611 template <typename TypeTag>
612 void
615 WellStateType& well_state,
616 DeferredLogger& deferred_logger)
617 {
618 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
619
620 // We assemble the well equations, then we check the convergence,
621 // which is why we do not put the assembleWellEq here.
622 try{
623 const BVectorWell dx_well = this->linSys_.solve();
624 updateWellState(simulator, dx_well, well_state, deferred_logger);
625 }
626 catch(const NumericalProblem& exp) {
627 // Add information about the well and log to deferred logger
628 // (Logging done inside of solve() method will only be seen if
629 // this is the process with rank zero)
630 deferred_logger.problem("In MultisegmentWell::solveEqAndUpdateWellState for well "
631 + this->name() +": "+exp.what());
632 throw;
633 }
634 }
635
636
637
638
639
640 template <typename TypeTag>
641 void
644 {
645 // We call this function on every process for the number_of_local_perforations_ on that process
646 // Each process updates the pressure for his perforations
647 for (int local_perf_index = 0; local_perf_index < this->number_of_local_perforations_; ++local_perf_index) {
648 // This variable loops over the number_of_local_perforations_ of *this* process, hence it is *local*.
649
650 std::vector<Scalar> kr(this->number_of_phases_, 0.0);
651 std::vector<Scalar> density(this->number_of_phases_, 0.0);
652
653 const int cell_idx = this->well_cells_[local_perf_index];
654 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
655 const auto& fs = intQuants.fluidState();
656
657 Scalar sum_kr = 0.;
658
659 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
660 const int water_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::waterPhaseIdx);
661 kr[water_pos] = intQuants.relativePermeability(FluidSystem::waterPhaseIdx).value();
662 sum_kr += kr[water_pos];
663 density[water_pos] = fs.density(FluidSystem::waterPhaseIdx).value();
664 }
665
666 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
667 const int oil_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::oilPhaseIdx);
668 kr[oil_pos] = intQuants.relativePermeability(FluidSystem::oilPhaseIdx).value();
669 sum_kr += kr[oil_pos];
670 density[oil_pos] = fs.density(FluidSystem::oilPhaseIdx).value();
671 }
672
673 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
674 const int gas_pos = FluidSystem::canonicalToActivePhaseIdx(FluidSystem::gasPhaseIdx);
675 kr[gas_pos] = intQuants.relativePermeability(FluidSystem::gasPhaseIdx).value();
676 sum_kr += kr[gas_pos];
677 density[gas_pos] = fs.density(FluidSystem::gasPhaseIdx).value();
678 }
679
680 assert(sum_kr != 0.);
681
682 // calculate the average density
683 Scalar average_density = 0.;
684 for (int p = 0; p < this->number_of_phases_; ++p) {
685 average_density += kr[p] * density[p];
686 }
687 average_density /= sum_kr;
688
689 this->cell_perforation_pressure_diffs_[local_perf_index] = this->gravity_ * average_density * this->cell_perforation_depth_diffs_[local_perf_index];
690 }
691 }
692
693
694
695
696
697 template <typename TypeTag>
698 void
701 DeferredLogger& deferred_logger)
702 {
703 for (int seg = 0; seg < this->numberOfSegments(); ++seg) {
704 // TODO: trying to reduce the times for the surfaceVolumeFraction calculation
705 const Scalar surface_volume = getSegmentSurfaceVolume(simulator, seg, deferred_logger).value();
706 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
707 segment_fluid_initial_[seg][comp_idx] = surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx).value();
708 }
709 }
710 }
711
712
713
714
715
716 template <typename TypeTag>
717 void
719 updateWellState(const Simulator& simulator,
720 const BVectorWell& dwells,
721 WellStateType& well_state,
722 DeferredLogger& deferred_logger,
723 const Scalar relaxation_factor)
724 {
725 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
726
727 const Scalar dFLimit = this->param_.dwell_fraction_max_;
728 const Scalar max_pressure_change = this->param_.max_pressure_change_ms_wells_;
729 const bool stop_or_zero_rate_target =
730 this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
731 this->primary_variables_.updateNewton(dwells,
732 relaxation_factor,
733 dFLimit,
734 stop_or_zero_rate_target,
735 max_pressure_change);
736
737 const auto& summary_state = simulator.vanguard().summaryState();
738 this->primary_variables_.copyToWellState(*this, getRefDensity(),
739 well_state,
740 summary_state,
741 deferred_logger);
742
743 {
744 auto& ws = well_state.well(this->index_of_well_);
745 this->segments_.copyPhaseDensities(ws.segments);
746 }
747
748 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.well(this->index_of_well_));
749 }
750
751
752
753
754
755 template <typename TypeTag>
756 void
759 const WellStateType& well_state,
760 DeferredLogger& deferred_logger)
761 {
762 updatePrimaryVariables(simulator, well_state, deferred_logger);
763 computePerfCellPressDiffs(simulator);
764 computeInitialSegmentFluids(simulator, deferred_logger);
765 }
766
767
768
769
770
771 template<typename TypeTag>
772 void
774 updateProductivityIndex(const Simulator& simulator,
775 const WellProdIndexCalculator<Scalar>& wellPICalc,
776 WellStateType& well_state,
777 DeferredLogger& deferred_logger) const
778 {
779 auto fluidState = [&simulator, this](const int local_perf_index)
780 {
781 const auto cell_idx = this->well_cells_[local_perf_index];
782 return simulator.model()
783 .intensiveQuantities(cell_idx, /*timeIdx=*/ 0).fluidState();
784 };
785
786 const int np = this->number_of_phases_;
787 auto setToZero = [np](Scalar* x) -> void
788 {
789 std::fill_n(x, np, 0.0);
790 };
791
792 auto addVector = [np](const Scalar* src, Scalar* dest) -> void
793 {
794 std::transform(src, src + np, dest, dest, std::plus<>{});
795 };
796
797 auto& ws = well_state.well(this->index_of_well_);
798 auto& perf_data = ws.perf_data;
799 auto* connPI = perf_data.prod_index.data();
800 auto* wellPI = ws.productivity_index.data();
801
802 setToZero(wellPI);
803
804 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
805 auto subsetPerfID = 0;
806
807 for ( const auto& perf : *this->perf_data_){
808 auto allPerfID = perf.ecl_index;
809
810 auto connPICalc = [&wellPICalc, allPerfID](const Scalar mobility) -> Scalar
811 {
812 return wellPICalc.connectionProdIndStandard(allPerfID, mobility);
813 };
814
815 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
816 // The subsetPerfID loops over 0 .. this->perf_data_->size().
817 // *(this->perf_data_) contains info about the local processes only,
818 // hence subsetPerfID is a local perf id and we can call getMobility
819 // as well as fluidState directly with that.
820 getMobility(simulator, static_cast<int>(subsetPerfID), mob, deferred_logger);
821
822 const auto& fs = fluidState(subsetPerfID);
823 setToZero(connPI);
824
825 if (this->isInjector()) {
826 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
827 mob, connPI, deferred_logger);
828 }
829 else { // Production or zero flow rate
830 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
831 }
832
833 addVector(connPI, wellPI);
834
835 ++subsetPerfID;
836 connPI += np;
837 }
838
839 // Sum with communication in case of distributed well.
840 const auto& comm = this->parallel_well_info_.communication();
841 if (comm.size() > 1) {
842 comm.sum(wellPI, np);
843 }
844
845 assert (static_cast<int>(subsetPerfID) == this->number_of_local_perforations_ &&
846 "Internal logic error in processing connections for PI/II");
847 }
848
849
850
851
852
853 template<typename TypeTag>
856 connectionDensity(const int globalConnIdx,
857 [[maybe_unused]] const int openConnIdx) const
858 {
859 // Simple approximation: Mixture density at reservoir connection is
860 // mixture density at connection's segment.
861
862 const auto segNum = this->wellEcl()
863 .getConnections()[globalConnIdx].segment();
864
865 const auto segIdx = this->wellEcl()
866 .getSegments().segmentNumberToIndex(segNum);
867
868 return this->segments_.density(segIdx).value();
869 }
870
871
872
873
874
875 template<typename TypeTag>
876 void
879 {
880 if (this->number_of_local_perforations_ == 0) {
881 // If there are no open perforations on this process, there are no contributions to the jacobian.
882 return;
883 }
884 this->linSys_.extract(jacobian);
885 }
886
887
888 template<typename TypeTag>
889 void
892 const BVector& weights,
893 const int pressureVarIndex,
894 const bool use_well_weights,
895 const WellStateType& well_state) const
896 {
897 if (this->number_of_local_perforations_ == 0) {
898 // If there are no open perforations on this process, there are no contributions the cpr pressure matrix.
899 return;
900 }
901 // Add the pressure contribution to the cpr system for the well
902 this->linSys_.extractCPRPressureMatrix(jacobian,
903 weights,
904 pressureVarIndex,
905 use_well_weights,
906 *this,
907 this->SPres,
908 well_state);
909 }
910
911
912 template<typename TypeTag>
913 template<class Value>
914 void
916 computePerfRate(const Value& pressure_cell,
917 const Value& rs,
918 const Value& rv,
919 const std::vector<Value>& b_perfcells,
920 const std::vector<Value>& mob_perfcells,
921 const std::vector<Value>& Tw,
922 const int perf,
923 const Value& segment_pressure,
924 const Value& segment_density,
925 const bool& allow_cf,
926 const std::vector<Value>& cmix_s,
927 std::vector<Value>& cq_s,
928 Value& perf_press,
929 PerforationRates<Scalar>& perf_rates,
930 DeferredLogger& deferred_logger) const
931 {
932 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
933 if (local_perf_index < 0) // then the perforation is not on this process
934 return;
935
936 // pressure difference between the segment and the perforation
937 const Value perf_seg_press_diff = this->gravity() * segment_density *
938 this->segments_.local_perforation_depth_diff(local_perf_index);
939 // pressure difference between the perforation and the grid cell
940 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
941
942 // perforation pressure is the wellbore pressure corrected to perforation depth
943 // (positive sign due to convention in segments_.local_perforation_depth_diff() )
944 perf_press = segment_pressure + perf_seg_press_diff;
945
946 // cell pressure corrected to perforation depth
947 const Value cell_press_at_perf = pressure_cell - cell_perf_press_diff;
948
949 // Pressure drawdown (also used to determine direction of flow)
950 const Value drawdown = cell_press_at_perf - perf_press;
951
952 // producing perforations
953 if (drawdown > 0.0) {
954 // Do nothing if crossflow is not allowed
955 if (!allow_cf && this->isInjector()) {
956 return;
957 }
958
959 // compute component volumetric rates at standard conditions
960 for (int comp_idx = 0; comp_idx < this->numConservationQuantities(); ++comp_idx) {
961 const Value cq_p = - Tw[comp_idx] * (mob_perfcells[comp_idx] * drawdown);
962 cq_s[comp_idx] = b_perfcells[comp_idx] * cq_p;
963 }
964
965 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
966 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
967 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
968 const Value cq_s_oil = cq_s[oilCompIdx];
969 const Value cq_s_gas = cq_s[gasCompIdx];
970 cq_s[gasCompIdx] += rs * cq_s_oil;
971 cq_s[oilCompIdx] += rv * cq_s_gas;
972 }
973 } else { // injecting perforations
974 // Do nothing if crossflow is not allowed
975 if (!allow_cf && this->isProducer()) {
976 return;
977 }
978
979 // for injecting perforations, we use total mobility
980 Value total_mob = mob_perfcells[0];
981 for (int comp_idx = 1; comp_idx < this->numConservationQuantities(); ++comp_idx) {
982 total_mob += mob_perfcells[comp_idx];
983 }
984
985 // compute volume ratio between connection and at standard conditions
986 Value volume_ratio = 0.0;
987 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
988 const unsigned waterCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
989 volume_ratio += cmix_s[waterCompIdx] / b_perfcells[waterCompIdx];
990 }
991
992 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
993 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
994 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
995
996 // Incorporate RS/RV factors if both oil and gas active
997 // TODO: not sure we use rs rv from the perforation cells when handling injecting perforations
998 // basically, for injecting perforations, the wellbore is the upstreaming side.
999 const Value d = 1.0 - rv * rs;
1000
1001 if (getValue(d) == 0.0) {
1002 OPM_DEFLOG_PROBLEM(NumericalProblem,
1003 fmt::format("Zero d value obtained for well {} "
1004 "during flux calculation with rs {} and rv {}",
1005 this->name(), rs, rv),
1006 deferred_logger);
1007 }
1008
1009 const Value tmp_oil = (cmix_s[oilCompIdx] - rv * cmix_s[gasCompIdx]) / d;
1010 volume_ratio += tmp_oil / b_perfcells[oilCompIdx];
1011
1012 const Value tmp_gas = (cmix_s[gasCompIdx] - rs * cmix_s[oilCompIdx]) / d;
1013 volume_ratio += tmp_gas / b_perfcells[gasCompIdx];
1014 } else { // not having gas and oil at the same time
1015 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
1016 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1017 volume_ratio += cmix_s[oilCompIdx] / b_perfcells[oilCompIdx];
1018 }
1019 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1020 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1021 volume_ratio += cmix_s[gasCompIdx] / b_perfcells[gasCompIdx];
1022 }
1023 }
1024 // injecting connections total volumerates at standard conditions
1025 for (int componentIdx = 0; componentIdx < this->numConservationQuantities(); ++componentIdx) {
1026 const Value cqt_i = - Tw[componentIdx] * (total_mob * drawdown);
1027 Value cqt_is = cqt_i / volume_ratio;
1028 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
1029 }
1030 } // end for injection perforations
1031
1032 // calculating the perforation solution gas rate and solution oil rates
1033 if (this->isProducer()) {
1034 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1035 const unsigned oilCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1036 const unsigned gasCompIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1037 // TODO: the formulations here remain to be tested with cases with strong crossflow through production wells
1038 // s means standard condition, r means reservoir condition
1039 // q_os = q_or * b_o + rv * q_gr * b_g
1040 // q_gs = q_gr * g_g + rs * q_or * b_o
1041 // d = 1.0 - rs * rv
1042 // q_or = 1 / (b_o * d) * (q_os - rv * q_gs)
1043 // q_gr = 1 / (b_g * d) * (q_gs - rs * q_os)
1044
1045 const Scalar d = 1.0 - getValue(rv) * getValue(rs);
1046 // vaporized oil into gas
1047 // rv * q_gr * b_g = rv * (q_gs - rs * q_os) / d
1048 perf_rates.vap_oil = getValue(rv) * (getValue(cq_s[gasCompIdx]) - getValue(rs) * getValue(cq_s[oilCompIdx])) / d;
1049 // dissolved of gas in oil
1050 // rs * q_or * b_o = rs * (q_os - rv * q_gs) / d
1051 perf_rates.dis_gas = getValue(rs) * (getValue(cq_s[oilCompIdx]) - getValue(rv) * getValue(cq_s[gasCompIdx])) / d;
1052 }
1053 }
1054 }
1055
1056 template <typename TypeTag>
1057 template<class Value>
1058 void
1060 computePerfRate(const IntensiveQuantities& int_quants,
1061 const std::vector<Value>& mob_perfcells,
1062 const std::vector<Value>& Tw,
1063 const int seg,
1064 const int perf,
1065 const Value& segment_pressure,
1066 const bool& allow_cf,
1067 std::vector<Value>& cq_s,
1068 Value& perf_press,
1069 PerforationRates<Scalar>& perf_rates,
1070 DeferredLogger& deferred_logger) const
1071
1072 {
1073 auto obtain = [this](const Eval& value)
1074 {
1075 if constexpr (std::is_same_v<Value, Scalar>) {
1076 static_cast<void>(this); // suppress clang warning
1077 return getValue(value);
1078 } else {
1079 return this->extendEval(value);
1080 }
1081 };
1082 auto obtainN = [](const auto& value)
1083 {
1084 if constexpr (std::is_same_v<Value, Scalar>) {
1085 return getValue(value);
1086 } else {
1087 return value;
1088 }
1089 };
1090 const auto& fs = int_quants.fluidState();
1091
1092 const Value pressure_cell = obtain(this->getPerfCellPressure(fs));
1093 const Value rs = obtain(fs.Rs());
1094 const Value rv = obtain(fs.Rv());
1095
1096 // not using number_of_phases_ because of solvent
1097 std::vector<Value> b_perfcells(this->num_conservation_quantities_, 0.0);
1098
1099 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
1100 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1101 continue;
1102 }
1103
1104 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1105 b_perfcells[compIdx] = obtain(fs.invB(phaseIdx));
1106 }
1107
1108 std::vector<Value> cmix_s(this->numConservationQuantities(), 0.0);
1109 for (int comp_idx = 0; comp_idx < this->numConservationQuantities(); ++comp_idx) {
1110 cmix_s[comp_idx] = obtainN(this->primary_variables_.surfaceVolumeFraction(seg, comp_idx));
1111 }
1112
1113 this->computePerfRate(pressure_cell,
1114 rs,
1115 rv,
1116 b_perfcells,
1117 mob_perfcells,
1118 Tw,
1119 perf,
1120 segment_pressure,
1121 obtainN(this->segments_.density(seg)),
1122 allow_cf,
1123 cmix_s,
1124 cq_s,
1125 perf_press,
1126 perf_rates,
1127 deferred_logger);
1128 }
1129
1130 template <typename TypeTag>
1131 void
1133 computeSegmentFluidProperties(const Simulator& simulator, DeferredLogger& deferred_logger)
1134 {
1135 // TODO: the concept of phases and components are rather confusing in this function.
1136 // needs to be addressed sooner or later.
1137
1138 // get the temperature for later use. It is only useful when we are not handling
1139 // thermal related simulation
1140 // basically, it is a single value for all the segments
1141
1142 EvalWell temperature;
1143 EvalWell saltConcentration;
1144 // not sure how to handle the pvt region related to segment
1145 // for the current approach, we use the pvt region of the first perforated cell
1146 // although there are some text indicating using the pvt region of the lowest
1147 // perforated cell
1148 // TODO: later to investigate how to handle the pvt region
1149
1150 auto info = this->getFirstPerforationFluidStateInfo(simulator);
1151 temperature.setValue(std::get<0>(info));
1152 saltConcentration = this->extendEval(std::get<1>(info));
1153
1154 this->segments_.computeFluidProperties(temperature,
1155 saltConcentration,
1156 this->primary_variables_,
1157 deferred_logger);
1158 }
1159
1160 template<typename TypeTag>
1161 template<class Value>
1162 void
1164 getTransMult(Value& trans_mult,
1165 const Simulator& simulator,
1166 const int cell_idx) const
1167 {
1168 auto obtain = [this](const Eval& value)
1169 {
1170 if constexpr (std::is_same_v<Value, Scalar>) {
1171 static_cast<void>(this); // suppress clang warning
1172 return getValue(value);
1173 } else {
1174 return this->extendEval(value);
1175 }
1176 };
1177 WellInterface<TypeTag>::getTransMult(trans_mult, simulator, cell_idx, obtain);
1178 }
1179
1180 template <typename TypeTag>
1181 template<class Value>
1182 void
1184 getMobility(const Simulator& simulator,
1185 const int local_perf_index,
1186 std::vector<Value>& mob,
1187 DeferredLogger& deferred_logger) const
1188 {
1189 auto obtain = [this](const Eval& value)
1190 {
1191 if constexpr (std::is_same_v<Value, Scalar>) {
1192 static_cast<void>(this); // suppress clang warning
1193 return getValue(value);
1194 } else {
1195 return this->extendEval(value);
1196 }
1197 };
1198
1199 WellInterface<TypeTag>::getMobility(simulator, local_perf_index, mob, obtain, deferred_logger);
1200
1201 if (this->isInjector() && this->well_ecl_.getInjMultMode() != Well::InjMultMode::NONE) {
1202 const auto perf_ecl_index = this->perforationData()[local_perf_index].ecl_index;
1203 const Connection& con = this->well_ecl_.getConnections()[perf_ecl_index];
1204 const int seg = this->segmentNumberToIndex(con.segment());
1205 // from the reference results, it looks like MSW uses segment pressure instead of BHP here
1206 // Note: this is against the documented definition.
1207 // we can change this depending on what we want
1208 const Scalar segment_pres = this->primary_variables_.getSegmentPressure(seg).value();
1209 const Scalar perf_seg_press_diff = this->gravity() * this->segments_.density(seg).value()
1210 * this->segments_.local_perforation_depth_diff(local_perf_index);
1211 const Scalar perf_press = segment_pres + perf_seg_press_diff;
1212 const Scalar multiplier = this->getInjMult(local_perf_index, segment_pres, perf_press, deferred_logger);
1213 for (std::size_t i = 0; i < mob.size(); ++i) {
1214 mob[i] *= multiplier;
1215 }
1216 }
1217 }
1218
1219
1220
1221 template<typename TypeTag>
1224 getRefDensity() const
1225 {
1226 return this->segments_.getRefDensity();
1227 }
1228
1229 template<typename TypeTag>
1230 void
1232 checkOperabilityUnderBHPLimit(const WellStateType& /*well_state*/,
1233 const Simulator& simulator,
1234 DeferredLogger& deferred_logger)
1235 {
1236 const auto& summaryState = simulator.vanguard().summaryState();
1237 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1238 // Crude but works: default is one atmosphere.
1239 // TODO: a better way to detect whether the BHP is defaulted or not
1240 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
1241 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
1242 // if the BHP limit is not defaulted or the well does not have a THP limit
1243 // we need to check the BHP limit
1244 Scalar total_ipr_mass_rate = 0.0;
1245 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
1246 {
1247 if (!FluidSystem::phaseIsActive(phaseIdx)) {
1248 continue;
1249 }
1250
1251 const unsigned compIdx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
1252 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
1253
1254 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
1255 total_ipr_mass_rate += ipr_rate * rho;
1256 }
1257 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
1258 this->operability_status_.operable_under_only_bhp_limit = false;
1259 }
1260
1261 // checking whether running under BHP limit will violate THP limit
1262 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
1263 // option 1: calculate well rates based on the BHP limit.
1264 // option 2: stick with the above IPR curve
1265 // we use IPR here
1266 std::vector<Scalar> well_rates_bhp_limit;
1267 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
1268
1269 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1270 const Scalar thp = WellBhpThpCalculator(*this).calculateThpFromBhp(well_rates_bhp_limit,
1271 bhp_limit,
1272 this->getRefDensity(),
1273 this->wellEcl().alq_value(summaryState),
1274 thp_limit,
1275 deferred_logger);
1276 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
1277 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
1278 }
1279 }
1280 } else {
1281 // defaulted BHP and there is a THP constraint
1282 // default BHP limit is about 1 atm.
1283 // when applied the hydrostatic pressure correction dp,
1284 // most likely we get a negative value (bhp + dp)to search in the VFP table,
1285 // which is not desirable.
1286 // we assume we can operate under defaulted BHP limit and will violate the THP limit
1287 // when operating under defaulted BHP limit.
1288 this->operability_status_.operable_under_only_bhp_limit = true;
1289 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
1290 }
1291 }
1292
1293
1294
1295 template<typename TypeTag>
1296 void
1298 updateIPR(const Simulator& simulator, DeferredLogger& deferred_logger) const
1299 {
1300 // TODO: not handling solvent related here for now
1301
1302 // initialize all the values to be zero to begin with
1303 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
1304 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
1305
1306 const int nseg = this->numberOfSegments();
1307 std::vector<Scalar> seg_dp(nseg, 0.0);
1308 for (int seg = 0; seg < nseg; ++seg) {
1309 // calculating the perforation rate for each perforation that belongs to this segment
1310 const Scalar dp = this->getSegmentDp(seg,
1311 this->segments_.density(seg).value(),
1312 seg_dp);
1313 seg_dp[seg] = dp;
1314 for (const int perf : this->segments_.perforations()[seg]) {
1315 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1316 if (local_perf_index < 0) // then the perforation is not on this process
1317 continue;
1318 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
1319
1320 // TODO: maybe we should store the mobility somewhere, so that we only need to calculate it one per iteration
1321 getMobility(simulator, local_perf_index, mob, deferred_logger);
1322
1323 const int cell_idx = this->well_cells_[local_perf_index];
1324 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1325 const auto& fs = int_quantities.fluidState();
1326 // pressure difference between the segment and the perforation
1327 const Scalar perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
1328 // pressure difference between the perforation and the grid cell
1329 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
1330 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
1331
1332 // calculating the b for the connection
1333 std::vector<Scalar> b_perf(this->num_conservation_quantities_);
1334 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
1335 if (!FluidSystem::phaseIsActive(phase)) {
1336 continue;
1337 }
1338 const unsigned comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phase));
1339 b_perf[comp_idx] = fs.invB(phase).value();
1340 }
1341
1342 // the pressure difference between the connection and BHP
1343 const Scalar h_perf = cell_perf_press_diff + perf_seg_press_diff + dp;
1344 const Scalar pressure_diff = pressure_cell - h_perf;
1345
1346 // do not take into consideration the crossflow here.
1347 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
1348 deferred_logger.debug("CROSSFLOW_IPR",
1349 "cross flow found when updateIPR for well " + this->name());
1350 }
1351
1352 // the well index associated with the connection
1353 Scalar trans_mult(0.0);
1354 getTransMult(trans_mult, simulator, cell_idx);
1355 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1356 std::vector<Scalar> tw_perf(this->num_conservation_quantities_, this->well_index_[perf] * trans_mult);
1357 this->getTw(tw_perf, local_perf_index, int_quantities, trans_mult, wellstate_nupcol);
1358 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
1359 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
1360 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1361 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
1362 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
1363 ipr_b_perf[comp_idx] += tw_mob;
1364 }
1365
1366 // we need to handle the rs and rv when both oil and gas are present
1367 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
1368 const unsigned oil_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::oilCompIdx);
1369 const unsigned gas_comp_idx = FluidSystem::canonicalToActiveCompIdx(FluidSystem::gasCompIdx);
1370 const Scalar rs = (fs.Rs()).value();
1371 const Scalar rv = (fs.Rv()).value();
1372
1373 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
1374 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
1375
1376 ipr_a_perf[gas_comp_idx] += dis_gas_a;
1377 ipr_a_perf[oil_comp_idx] += vap_oil_a;
1378
1379 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
1380 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
1381
1382 ipr_b_perf[gas_comp_idx] += dis_gas_b;
1383 ipr_b_perf[oil_comp_idx] += vap_oil_b;
1384 }
1385
1386 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
1387 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
1388 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
1389 }
1390 }
1391 }
1392 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
1393 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
1394 }
1395
1396 template<typename TypeTag>
1397 void
1399 updateIPRImplicit(const Simulator& simulator,
1400 WellStateType& well_state,
1401 DeferredLogger& deferred_logger)
1402 {
1403 // Compute IPR based on *converged* well-equation:
1404 // For a component rate r the derivative dr/dbhp is obtained by
1405 // dr/dbhp = - (partial r/partial x) * inv(partial Eq/partial x) * (partial Eq/partial bhp_target)
1406 // where Eq(x)=0 is the well equation setup with bhp control and primary variables x
1407
1408 // We shouldn't have zero rates at this stage, but check
1409 bool zero_rates;
1410 auto rates = well_state.well(this->index_of_well_).surface_rates;
1411 zero_rates = true;
1412 for (std::size_t p = 0; p < rates.size(); ++p) {
1413 zero_rates &= rates[p] == 0.0;
1414 }
1415 auto& ws = well_state.well(this->index_of_well_);
1416 if (zero_rates) {
1417 const auto msg = fmt::format("updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
1418 deferred_logger.debug(msg);
1419 /*
1420 // could revert to standard approach here:
1421 updateIPR(simulator, deferred_logger);
1422 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
1423 const int idx = this->activeCompToActivePhaseIdx(comp_idx);
1424 ws.implicit_ipr_a[idx] = this->ipr_a_[comp_idx];
1425 ws.implicit_ipr_b[idx] = this->ipr_b_[comp_idx];
1426 }
1427 return;
1428 */
1429 }
1430 const auto& group_state = simulator.problem().wellModel().groupState();
1431
1432 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
1433 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
1434 //WellState well_state_copy = well_state;
1435 auto inj_controls = Well::InjectionControls(0);
1436 auto prod_controls = Well::ProductionControls(0);
1437 prod_controls.addControl(Well::ProducerCMode::BHP);
1438 prod_controls.bhp_limit = well_state.well(this->index_of_well_).bhp;
1439
1440 // Set current control to bhp, and bhp value in state, modify bhp limit in control object.
1441 const auto cmode = ws.production_cmode;
1442 ws.production_cmode = Well::ProducerCMode::BHP;
1443 const double dt = simulator.timeStepSize();
1444 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
1445
1446 BVectorWell rhs(this->numberOfSegments());
1447 rhs = 0.0;
1448 rhs[0][SPres] = -1.0;
1449
1450 const BVectorWell x_well = this->linSys_.solve(rhs);
1451 constexpr int num_eq = MSWEval::numWellEq;
1452 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx){
1453 const EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
1454 const int idx = FluidSystem::activeCompToActivePhaseIdx(comp_idx);
1455 for (size_t pvIdx = 0; pvIdx < num_eq; ++pvIdx) {
1456 // well primary variable derivatives in EvalWell start at position Indices::numEq
1457 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
1458 }
1459 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
1460 }
1461 // reset cmode
1462 ws.production_cmode = cmode;
1463 }
1464
1465 template<typename TypeTag>
1466 void
1469 const WellStateType& well_state,
1470 const WellGroupHelperType& wgHelper,
1471 DeferredLogger& deferred_logger)
1472 {
1473 const auto& summaryState = simulator.vanguard().summaryState();
1474 const auto obtain_bhp = this->isProducer()
1475 ? computeBhpAtThpLimitProd(
1476 well_state, simulator, wgHelper, summaryState, deferred_logger)
1477 : computeBhpAtThpLimitInj(simulator, wgHelper, summaryState, deferred_logger);
1478
1479 if (obtain_bhp) {
1480 this->operability_status_.can_obtain_bhp_with_thp_limit = true;
1481
1482 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1483 this->operability_status_.obey_bhp_limit_with_thp_limit = (*obtain_bhp >= bhp_limit);
1484
1485 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1486 if (this->isProducer() && *obtain_bhp < thp_limit) {
1487 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1488 + " bars is SMALLER than thp limit "
1489 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1490 + " bars as a producer for well " + this->name();
1491 deferred_logger.debug(msg);
1492 }
1493 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1494 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1495 + " bars is LARGER than thp limit "
1496 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1497 + " bars as a injector for well " + this->name();
1498 deferred_logger.debug(msg);
1499 }
1500 } else {
1501 // Shutting wells that can not find bhp value from thp
1502 // when under THP control
1503 this->operability_status_.can_obtain_bhp_with_thp_limit = false;
1504 this->operability_status_.obey_bhp_limit_with_thp_limit = false;
1505 if (!this->wellIsStopped()) {
1506 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1507 deferred_logger.debug(" could not find bhp value at thp limit "
1508 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1509 + " bar for well " + this->name() + ", the well might need to be closed ");
1510 }
1511 }
1512 }
1513
1514
1515
1516
1517
1518 template<typename TypeTag>
1519 bool
1521 iterateWellEqWithControl(const Simulator& simulator,
1522 const double dt,
1523 const Well::InjectionControls& inj_controls,
1524 const Well::ProductionControls& prod_controls,
1525 const WellGroupHelperType& wgHelper,
1526 WellStateType& well_state,
1527 DeferredLogger& deferred_logger)
1528 {
1529 const auto& group_state = wgHelper.groupState();
1530 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return true;
1531
1532 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1533
1534 {
1535 // getWellFiniteResiduals returns false for nan/inf residuals
1536 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1537 if(!isFinite)
1538 return false;
1539 }
1540
1541 updatePrimaryVariables(simulator, well_state, deferred_logger);
1542
1543 std::vector<std::vector<Scalar> > residual_history;
1544 std::vector<Scalar> measure_history;
1545 int it = 0;
1546 // relaxation factor
1547 Scalar relaxation_factor = 1.;
1548 bool converged = false;
1549 bool relax_convergence = false;
1550 this->regularize_ = false;
1551 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1552
1553 if (it > this->param_.strict_inner_iter_wells_) {
1554 relax_convergence = true;
1555 this->regularize_ = true;
1556 }
1557
1558 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls,
1559 well_state, group_state, deferred_logger);
1560
1561 const auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
1562 if (report.converged()) {
1563 converged = true;
1564 break;
1565 }
1566
1567 {
1568 // getFinteWellResiduals returns false for nan/inf residuals
1569 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1570 if (!isFinite)
1571 return false;
1572
1573 residual_history.push_back(residuals);
1574 measure_history.push_back(this->getResidualMeasureValue(well_state,
1575 residual_history[it],
1576 this->param_.tolerance_wells_,
1577 this->param_.tolerance_pressure_ms_wells_,
1578 deferred_logger) );
1579 }
1580 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1581 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1582 // try last attempt with relaxed tolerances
1583 const auto reportStag = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, true);
1584 if (reportStag.converged()) {
1585 converged = true;
1586 std::string message = fmt::format("Well stagnates/oscillates but {} manages to get converged with relaxed tolerances in {} inner iterations."
1587 ,this->name(), it);
1588 deferred_logger.debug(message);
1589 } else {
1590 converged = false;
1591 }
1592 break;
1593 }
1594
1595 BVectorWell dx_well;
1596 try{
1597 dx_well = this->linSys_.solve();
1598 updateWellState(simulator, dx_well, well_state, deferred_logger, relaxation_factor);
1599 }
1600 catch(const NumericalProblem& exp) {
1601 // Add information about the well and log to deferred logger
1602 // (Logging done inside of solve() method will only be seen if
1603 // this is the process with rank zero)
1604 deferred_logger.problem("In MultisegmentWell::iterateWellEqWithControl for well "
1605 + this->name() +": "+exp.what());
1606 throw;
1607 }
1608 }
1609
1610 // TODO: we should decide whether to keep the updated well_state, or recover to use the old well_state
1611 if (converged) {
1612 std::ostringstream sstr;
1613 sstr << " Well " << this->name() << " converged in " << it << " inner iterations.";
1614 if (relax_convergence)
1615 sstr << " (A relaxed tolerance was used after "<< this->param_.strict_inner_iter_wells_ << " iterations)";
1616
1617 // Output "converged in 0 inner iterations" messages only at
1618 // elevated verbosity levels.
1619 deferred_logger.debug(sstr.str(), OpmLog::defaultDebugVerbosityLevel + (it == 0));
1620 } else {
1621 std::ostringstream sstr;
1622 sstr << " Well " << this->name() << " did not converge in " << it << " inner iterations.";
1623#define EXTRA_DEBUG_MSW 0
1624#if EXTRA_DEBUG_MSW
1625 sstr << "***** Outputting the residual history for well " << this->name() << " during inner iterations:";
1626 for (int i = 0; i < it; ++i) {
1627 const auto& residual = residual_history[i];
1628 sstr << " residual at " << i << "th iteration ";
1629 for (const auto& res : residual) {
1630 sstr << " " << res;
1631 }
1632 sstr << " " << measure_history[i] << " \n";
1633 }
1634#endif
1635#undef EXTRA_DEBUG_MSW
1636 deferred_logger.debug(sstr.str());
1637 }
1638
1639 return converged;
1640 }
1641
1642
1643 template<typename TypeTag>
1644 bool
1646 iterateWellEqWithSwitching(const Simulator& simulator,
1647 const double dt,
1648 const Well::InjectionControls& inj_controls,
1649 const Well::ProductionControls& prod_controls,
1650 const WellGroupHelperType& wgHelper,
1651 WellStateType& well_state,
1652 DeferredLogger& deferred_logger,
1653 const bool fixed_control /*false*/,
1654 const bool fixed_status /*false*/)
1655 {
1656 const auto& group_state = wgHelper.groupState();
1657 const int max_iter_number = this->param_.max_inner_iter_ms_wells_;
1658
1659 {
1660 // getWellFiniteResiduals returns false for nan/inf residuals
1661 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1662 if(!isFinite)
1663 return false;
1664 }
1665
1666 updatePrimaryVariables(simulator, well_state, deferred_logger);
1667
1668 std::vector<std::vector<Scalar> > residual_history;
1669 std::vector<Scalar> measure_history;
1670 int it = 0;
1671 // relaxation factor
1672 Scalar relaxation_factor = 1.;
1673 bool converged = false;
1674 bool relax_convergence = false;
1675 this->regularize_ = false;
1676 const auto& summary_state = simulator.vanguard().summaryState();
1677
1678 // Always take a few (more than one) iterations after a switch before allowing a new switch
1679 // The optimal number here is subject to further investigation, but it has been observerved
1680 // that unless this number is >1, we may get stuck in a cycle
1681 const int min_its_after_switch = 3;
1682 // We also want to restrict the number of status switches to avoid oscillation between STOP<->OPEN
1683 const int max_status_switch = this->param_.max_well_status_switch_inner_iter_;
1684 int its_since_last_switch = min_its_after_switch;
1685 int switch_count= 0;
1686 int status_switch_count = 0;
1687 // if we fail to solve eqs, we reset status/operability before leaving
1688 const auto well_status_orig = this->wellStatus_;
1689 const auto operability_orig = this->operability_status_;
1690 auto well_status_cur = well_status_orig;
1691 // don't allow opening wells that has a stopped well status
1692 const bool allow_open = well_state.well(this->index_of_well_).status == WellStatus::OPEN;
1693 // don't allow switcing for wells under zero rate target or requested fixed status and control
1694 const bool allow_switching = !this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) &&
1695 (!fixed_control || !fixed_status) && allow_open;
1696 bool final_check = false;
1697 // well needs to be set operable or else solving/updating of re-opened wells is skipped
1698 this->operability_status_.resetOperability();
1699 this->operability_status_.solvable = true;
1700
1701 for (; it < max_iter_number; ++it, ++debug_cost_counter_) {
1702 ++its_since_last_switch;
1703 if (allow_switching && its_since_last_switch >= min_its_after_switch && status_switch_count < max_status_switch){
1704 const Scalar wqTotal = this->primary_variables_.getWQTotal().value();
1705 bool changed = this->updateWellControlAndStatusLocalIteration(
1706 simulator, wgHelper, inj_controls, prod_controls, wqTotal,
1707 well_state, deferred_logger, fixed_control, fixed_status
1708 );
1709 if (changed) {
1710 its_since_last_switch = 0;
1711 ++switch_count;
1712 if (well_status_cur != this->wellStatus_) {
1713 well_status_cur = this->wellStatus_;
1714 status_switch_count++;
1715 }
1716 }
1717 if (!changed && final_check) {
1718 break;
1719 } else {
1720 final_check = false;
1721 }
1722 if (status_switch_count == max_status_switch) {
1723 this->wellStatus_ = well_status_orig;
1724 }
1725 }
1726
1727 if (it > this->param_.strict_inner_iter_wells_) {
1728 relax_convergence = true;
1729 this->regularize_ = true;
1730 }
1731
1732 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls,
1733 well_state, group_state, deferred_logger);
1734
1735
1736 const auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
1737 converged = report.converged();
1738 if (this->parallel_well_info_.communication().size() > 1 &&
1739 this->parallel_well_info_.communication().max(converged) != this->parallel_well_info_.communication().min(converged)) {
1740 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()));
1741 }
1742 if (converged) {
1743 // if equations are sufficiently linear they might converge in less than min_its_after_switch
1744 // in this case, make sure all constraints are satisfied before returning
1745 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
1746 final_check = true;
1747 its_since_last_switch = min_its_after_switch;
1748 } else {
1749 break;
1750 }
1751 }
1752
1753 // getFinteWellResiduals returns false for nan/inf residuals
1754 {
1755 const auto& [isFinite, residuals] = this->getFiniteWellResiduals(Base::B_avg_, deferred_logger);
1756 if (!isFinite) {
1757 converged = false; // Jump out of loop instead of returning to ensure operability status is recovered
1758 break;
1759 }
1760
1761 residual_history.push_back(residuals);
1762 }
1763
1764 if (!converged) {
1765 measure_history.push_back(this->getResidualMeasureValue(well_state,
1766 residual_history[it],
1767 this->param_.tolerance_wells_,
1768 this->param_.tolerance_pressure_ms_wells_,
1769 deferred_logger));
1770 bool min_relaxation_reached = this->update_relaxation_factor(measure_history, relaxation_factor, this->regularize_, deferred_logger);
1771 if (min_relaxation_reached || this->repeatedStagnation(measure_history, this->regularize_, deferred_logger)) {
1772 converged = false;
1773 break;
1774 }
1775 }
1776 try{
1777 const BVectorWell dx_well = this->linSys_.solve();
1778 updateWellState(simulator, dx_well, well_state, deferred_logger, relaxation_factor);
1779 }
1780 catch(const NumericalProblem& exp) {
1781 // Add information about the well and log to deferred logger
1782 // (Logging done inside of solve() method will only be seen if
1783 // this is the process with rank zero)
1784 deferred_logger.problem("In MultisegmentWell::iterateWellEqWithSwitching for well "
1785 + this->name() +": "+exp.what());
1786 throw;
1787 }
1788 }
1789
1790 if (converged) {
1791 if (allow_switching){
1792 // update operability if status change
1793 const bool is_stopped = this->wellIsStopped();
1794 if (this->wellHasTHPConstraints(summary_state)){
1795 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
1796 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
1797 } else {
1798 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
1799 }
1800 }
1801 std::string message = fmt::format(" Well {} converged in {} inner iterations ("
1802 "{} control/status switches).", this->name(), it, switch_count);
1803 if (relax_convergence) {
1804 message.append(fmt::format(" (A relaxed tolerance was used after {} iterations)",
1805 this->param_.strict_inner_iter_wells_));
1806 }
1807 deferred_logger.debug(message, OpmLog::defaultDebugVerbosityLevel + ((it == 0) && (switch_count == 0)));
1808 } else {
1809 this->wellStatus_ = well_status_orig;
1810 this->operability_status_ = operability_orig;
1811 const std::string message = fmt::format(" Well {} did not converge in {} inner iterations ("
1812 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
1813 deferred_logger.debug(message);
1814 this->primary_variables_.outputLowLimitPressureSegments(deferred_logger);
1815 }
1816
1817 return converged;
1818 }
1819
1820
1821 template<typename TypeTag>
1822 void
1825 const double dt,
1826 const Well::InjectionControls& inj_controls,
1827 const Well::ProductionControls& prod_controls,
1828 WellStateType& well_state,
1829 const GroupState<Scalar>& group_state,
1830 DeferredLogger& deferred_logger)
1831 {
1832 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1833
1834 // update the upwinding segments
1835 this->segments_.updateUpwindingSegments(this->primary_variables_);
1836
1837 // calculate the fluid properties needed.
1838 computeSegmentFluidProperties(simulator, deferred_logger);
1839
1840 // clear all entries
1841 this->linSys_.clear();
1842
1843 auto& ws = well_state.well(this->index_of_well_);
1844 ws.phase_mixing_rates.fill(0.0);
1845
1846 // for the black oil cases, there will be four equations,
1847 // the first three of them are the mass balance equations, the last one is the pressure equations.
1848 //
1849 // but for the top segment, the pressure equation will be the well control equation, and the other three will be the same.
1850
1851 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1852
1853 const int nseg = this->numberOfSegments();
1854
1855 const Scalar rhow = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1856 FluidSystem::referenceDensity( FluidSystem::waterPhaseIdx, Base::pvtRegionIdx() ) : 0.0;
1857 const unsigned watCompIdx = FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) ?
1858 FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx) : 0;
1859
1860 for (int seg = 0; seg < nseg; ++seg) {
1861 // calculating the perforation rate for each perforation that belongs to this segment
1862 const EvalWell seg_pressure = this->primary_variables_.getSegmentPressure(seg);
1863 auto& perf_data = ws.perf_data;
1864 auto& perf_rates = perf_data.phase_rates;
1865 auto& perf_press_state = perf_data.pressure;
1866 for (const int perf : this->segments_.perforations()[seg]) {
1867 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
1868 if (local_perf_index < 0) // then the perforation is not on this process
1869 continue;
1870 const int cell_idx = this->well_cells_[local_perf_index];
1871 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1872 std::vector<EvalWell> mob(this->num_conservation_quantities_, 0.0);
1873 getMobility(simulator, local_perf_index, mob, deferred_logger);
1874 EvalWell trans_mult(0.0);
1875 getTransMult(trans_mult, simulator, cell_idx);
1876 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1877 std::vector<EvalWell> Tw(this->num_conservation_quantities_, this->well_index_[local_perf_index] * trans_mult);
1878 this->getTw(Tw, local_perf_index, int_quants, trans_mult, wellstate_nupcol);
1879 std::vector<EvalWell> cq_s(this->num_conservation_quantities_, 0.0);
1880 EvalWell perf_press;
1881 PerforationRates<Scalar> perfRates;
1882 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
1883 allow_cf, cq_s, perf_press, perfRates, deferred_logger);
1884
1885 // updating the solution gas rate and solution oil rate
1886 if (this->isProducer()) {
1887 ws.phase_mixing_rates[ws.dissolved_gas] += perfRates.dis_gas;
1888 ws.phase_mixing_rates[ws.vaporized_oil] += perfRates.vap_oil;
1889 perf_data.phase_mixing_rates[local_perf_index][ws.dissolved_gas] = perfRates.dis_gas;
1890 perf_data.phase_mixing_rates[local_perf_index][ws.vaporized_oil] = perfRates.vap_oil;
1891 }
1892
1893 // store the perf pressure and rates
1894 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1895 perf_rates[local_perf_index*this->number_of_phases_ + FluidSystem::activeCompToActivePhaseIdx(comp_idx)] = cq_s[comp_idx].value();
1896 }
1897 perf_press_state[local_perf_index] = perf_press.value();
1898
1899 // mass rates, for now only water
1900 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
1901 perf_data.wat_mass_rates[local_perf_index] = cq_s[watCompIdx].value() * rhow;
1902 }
1903
1904 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1905 // the cq_s entering mass balance equations need to consider the efficiency factors.
1906 const EvalWell cq_s_effective = cq_s[comp_idx] * this->well_efficiency_factor_;
1907
1908 this->connectionRates_[local_perf_index][comp_idx] = Base::restrictEval(cq_s_effective);
1909
1911 assemblePerforationEq(seg, local_perf_index, comp_idx, cq_s_effective, this->linSys_);
1912 }
1913 }
1914 }
1915 // Accumulate dissolved gas and vaporized oil flow rates across all ranks sharing this well.
1916 {
1917 const auto& comm = this->parallel_well_info_.communication();
1918 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
1919 }
1920
1921 if (this->parallel_well_info_.communication().size() > 1) {
1922 // accumulate resWell_ and duneD_ in parallel to get effects of all perforations (might be distributed)
1923 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
1924 }
1925 for (int seg = 0; seg < nseg; ++seg) {
1926 // calculating the accumulation term
1927 // TODO: without considering the efficiency factor for now
1928 {
1929 const EvalWell segment_surface_volume = getSegmentSurfaceVolume(simulator, seg, deferred_logger);
1930
1931 // Add a regularization_factor to increase the accumulation term
1932 // This will make the system less stiff and help convergence for
1933 // difficult cases
1934 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
1935 // for each component
1936 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1937 const EvalWell accumulation_term = regularization_factor * (segment_surface_volume * this->primary_variables_.surfaceVolumeFraction(seg, comp_idx)
1938 - segment_fluid_initial_[seg][comp_idx]) / dt;
1940 assembleAccumulationTerm(seg, comp_idx, accumulation_term, this->linSys_);
1941 }
1942 }
1943 // considering the contributions due to flowing out from the segment
1944 {
1945 const int seg_upwind = this->segments_.upwinding_segment(seg);
1946 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1947 const EvalWell segment_rate =
1948 this->primary_variables_.getSegmentRateUpwinding(seg,
1949 seg_upwind,
1950 comp_idx) *
1951 this->well_efficiency_factor_;
1953 assembleOutflowTerm(seg, seg_upwind, comp_idx, segment_rate, this->linSys_);
1954 }
1955 }
1956
1957 // considering the contributions from the inlet segments
1958 {
1959 for (const int inlet : this->segments_.inlets()[seg]) {
1960 const int inlet_upwind = this->segments_.upwinding_segment(inlet);
1961 for (int comp_idx = 0; comp_idx < this->num_conservation_quantities_; ++comp_idx) {
1962 const EvalWell inlet_rate =
1963 this->primary_variables_.getSegmentRateUpwinding(inlet,
1964 inlet_upwind,
1965 comp_idx) *
1966 this->well_efficiency_factor_;
1968 assembleInflowTerm(seg, inlet, inlet_upwind, comp_idx, inlet_rate, this->linSys_);
1969 }
1970 }
1971 }
1972
1973 // the fourth equation, the pressure drop equation
1974 if (seg == 0) { // top segment, pressure equation is the control equation
1975 const auto& summaryState = simulator.vanguard().summaryState();
1976 const Schedule& schedule = simulator.vanguard().schedule();
1977 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1979 assembleControlEq(well_state,
1980 group_state,
1981 schedule,
1982 summaryState,
1983 inj_controls,
1984 prod_controls,
1985 getRefDensity(),
1986 this->primary_variables_,
1987 this->linSys_,
1988 stopped_or_zero_target,
1989 deferred_logger);
1990 } else {
1991 const UnitSystem& unit_system = simulator.vanguard().eclState().getDeckUnitSystem();
1992 const auto& summary_state = simulator.vanguard().summaryState();
1993 this->assemblePressureEq(seg, unit_system, well_state, summary_state, this->param_.use_average_density_ms_wells_, deferred_logger);
1994 }
1995 }
1996
1997 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
1998 this->linSys_.createSolver();
1999 }
2000
2001
2002
2003
2004 template<typename TypeTag>
2005 bool
2007 openCrossFlowAvoidSingularity(const Simulator& simulator) const
2008 {
2009 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
2010 }
2011
2012
2013 template<typename TypeTag>
2014 bool
2016 allDrawDownWrongDirection(const Simulator& simulator) const
2017 {
2018 bool all_drawdown_wrong_direction = true;
2019 const int nseg = this->numberOfSegments();
2020
2021 for (int seg = 0; seg < nseg; ++seg) {
2022 const EvalWell segment_pressure = this->primary_variables_.getSegmentPressure(seg);
2023 for (const int perf : this->segments_.perforations()[seg]) {
2024 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2025 if (local_perf_index < 0) // then the perforation is not on this process
2026 continue;
2027
2028 const int cell_idx = this->well_cells_[local_perf_index];
2029 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2030 const auto& fs = intQuants.fluidState();
2031
2032 // pressure difference between the segment and the perforation
2033 const EvalWell perf_seg_press_diff = this->segments_.getPressureDiffSegLocalPerf(seg, local_perf_index);
2034 // pressure difference between the perforation and the grid cell
2035 const Scalar cell_perf_press_diff = this->cell_perforation_pressure_diffs_[local_perf_index];
2036
2037 const Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2038 const Scalar perf_press = pressure_cell - cell_perf_press_diff;
2039 // Pressure drawdown (also used to determine direction of flow)
2040 // TODO: not 100% sure about the sign of the seg_perf_press_diff
2041 const EvalWell drawdown = perf_press - (segment_pressure + perf_seg_press_diff);
2042
2043 // for now, if there is one perforation can produce/inject in the correct
2044 // direction, we consider this well can still produce/inject.
2045 // TODO: it can be more complicated than this to cause wrong-signed rates
2046 if ( (drawdown < 0. && this->isInjector()) ||
2047 (drawdown > 0. && this->isProducer()) ) {
2048 all_drawdown_wrong_direction = false;
2049 break;
2050 }
2051 }
2052 }
2053 const auto& comm = this->parallel_well_info_.communication();
2054 if (comm.size() > 1)
2055 {
2056 all_drawdown_wrong_direction =
2057 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
2058 }
2059
2060 return all_drawdown_wrong_direction;
2061 }
2062
2063
2064
2065
2066 template<typename TypeTag>
2067 void
2069 updateWaterThroughput(const double /*dt*/, WellStateType& /*well_state*/) const
2070 {
2071 }
2072
2073
2074
2075
2076
2077 template<typename TypeTag>
2080 getSegmentSurfaceVolume(const Simulator& simulator,
2081 const int seg_idx,
2082 DeferredLogger& deferred_logger) const
2083 {
2084 EvalWell temperature;
2085 EvalWell saltConcentration;
2086
2087 auto info = this->getFirstPerforationFluidStateInfo(simulator);
2088 temperature.setValue(std::get<0>(info));
2089 saltConcentration = this->extendEval(std::get<1>(info));
2090
2091 return this->segments_.getSurfaceVolume(temperature,
2092 saltConcentration,
2093 this->primary_variables_,
2094 seg_idx,
2095 deferred_logger);
2096 }
2097
2098
2099 template<typename TypeTag>
2100 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2103 const Simulator& simulator,
2104 const WellGroupHelperType& wgHelper,
2105 const SummaryState& summary_state,
2106 DeferredLogger& deferred_logger) const
2107 {
2109 simulator,
2110 wgHelper,
2111 summary_state,
2112 this->getALQ(well_state),
2113 deferred_logger,
2114 /*iterate_if_no_solution */ true);
2115 }
2116
2117
2118
2119 template<typename TypeTag>
2120 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2123 const WellGroupHelperType& wgHelper,
2124 const SummaryState& summary_state,
2125 const Scalar alq_value,
2126 DeferredLogger& deferred_logger,
2127 bool iterate_if_no_solution) const
2128 {
2129 OPM_TIMEFUNCTION();
2130 // Make the frates() function.
2131 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2132 // Not solving the well equations here, which means we are
2133 // calculating at the current Fg/Fw values of the
2134 // well. This does not matter unless the well is
2135 // crossflowing, and then it is likely still a good
2136 // approximation.
2137 std::vector<Scalar> rates(3);
2138 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2139 return rates;
2140 };
2141
2142 auto bhpAtLimit = WellBhpThpCalculator(*this).
2143 computeBhpAtThpLimitProd(frates,
2144 summary_state,
2145 this->maxPerfPress(simulator),
2146 this->getRefDensity(),
2147 alq_value,
2148 this->getTHPConstraint(summary_state),
2149 deferred_logger);
2150
2151 if (bhpAtLimit)
2152 return bhpAtLimit;
2153
2154 if (!iterate_if_no_solution)
2155 return std::nullopt;
2156
2157 auto fratesIter = [this, &simulator, &wgHelper, &deferred_logger](const Scalar bhp) {
2158 // Solver the well iterations to see if we are
2159 // able to get a solution with an update
2160 // solution
2161 std::vector<Scalar> rates(3);
2162 computeWellRatesWithBhpIterations(simulator, bhp, wgHelper, rates, deferred_logger);
2163 return rates;
2164 };
2165
2166 return WellBhpThpCalculator(*this).
2167 computeBhpAtThpLimitProd(fratesIter,
2168 summary_state,
2169 this->maxPerfPress(simulator),
2170 this->getRefDensity(),
2171 alq_value,
2172 this->getTHPConstraint(summary_state),
2173 deferred_logger);
2174 }
2175
2176 template<typename TypeTag>
2177 std::optional<typename MultisegmentWell<TypeTag>::Scalar>
2179 computeBhpAtThpLimitInj(const Simulator& simulator,
2180 const WellGroupHelperType& wgHelper,
2181 const SummaryState& summary_state,
2182 DeferredLogger& deferred_logger) const
2183 {
2184 // Make the frates() function.
2185 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2186 // Not solving the well equations here, which means we are
2187 // calculating at the current Fg/Fw values of the
2188 // well. This does not matter unless the well is
2189 // crossflowing, and then it is likely still a good
2190 // approximation.
2191 std::vector<Scalar> rates(3);
2192 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2193 return rates;
2194 };
2195
2196 auto bhpAtLimit = WellBhpThpCalculator(*this).
2197 computeBhpAtThpLimitInj(frates,
2198 summary_state,
2199 this->getRefDensity(),
2200 0.05,
2201 100,
2202 false,
2203 deferred_logger);
2204
2205 if (bhpAtLimit)
2206 return bhpAtLimit;
2207
2208 auto fratesIter = [this, &simulator, &wgHelper, &deferred_logger](const Scalar bhp) {
2209 // Solver the well iterations to see if we are
2210 // able to get a solution with an update
2211 // solution
2212 std::vector<Scalar> rates(3);
2213 computeWellRatesWithBhpIterations(simulator, bhp, wgHelper, rates, deferred_logger);
2214 return rates;
2215 };
2216
2217 return WellBhpThpCalculator(*this).
2218 computeBhpAtThpLimitInj(fratesIter,
2219 summary_state,
2220 this->getRefDensity(),
2221 0.05,
2222 100,
2223 false,
2224 deferred_logger);
2225 }
2226
2227
2228
2229
2230
2231 template<typename TypeTag>
2234 maxPerfPress(const Simulator& simulator) const
2235 {
2236 Scalar max_pressure = 0.0;
2237 const int nseg = this->numberOfSegments();
2238 for (int seg = 0; seg < nseg; ++seg) {
2239 for (const int perf : this->segments_.perforations()[seg]) {
2240 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2241 if (local_perf_index < 0) // then the perforation is not on this process
2242 continue;
2243
2244 const int cell_idx = this->well_cells_[local_perf_index];
2245 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2246 const auto& fs = int_quants.fluidState();
2247 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2248 max_pressure = std::max(max_pressure, pressure_cell);
2249 }
2250 }
2251 max_pressure = this->parallel_well_info_.communication().max(max_pressure);
2252 return max_pressure;
2253 }
2254
2255
2256
2257
2258
2259 template<typename TypeTag>
2260 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2262 computeCurrentWellRates(const Simulator& simulator,
2263 DeferredLogger& deferred_logger) const
2264 {
2265 // Calculate the rates that follow from the current primary variables.
2266 std::vector<Scalar> well_q_s(this->num_conservation_quantities_, 0.0);
2267 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2268 const int nseg = this->numberOfSegments();
2269 for (int seg = 0; seg < nseg; ++seg) {
2270 // calculating the perforation rate for each perforation that belongs to this segment
2271 const Scalar seg_pressure = getValue(this->primary_variables_.getSegmentPressure(seg));
2272 for (const int perf : this->segments_.perforations()[seg]) {
2273 const int local_perf_index = this->pw_info_.activePerfToLocalPerf(perf);
2274 if (local_perf_index < 0) // then the perforation is not on this process
2275 continue;
2276
2277 const int cell_idx = this->well_cells_[local_perf_index];
2278 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2279 std::vector<Scalar> mob(this->num_conservation_quantities_, 0.0);
2280 getMobility(simulator, local_perf_index, mob, deferred_logger);
2281 Scalar trans_mult(0.0);
2282 getTransMult(trans_mult, simulator, cell_idx);
2283 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2284 std::vector<Scalar> Tw(this->num_conservation_quantities_, this->well_index_[local_perf_index] * trans_mult);
2285 this->getTw(Tw, local_perf_index, int_quants, trans_mult, wellstate_nupcol);
2286 std::vector<Scalar> cq_s(this->num_conservation_quantities_, 0.0);
2287 Scalar perf_press = 0.0;
2288 PerforationRates<Scalar> perf_rates;
2289 computePerfRate(int_quants, mob, Tw, seg, perf, seg_pressure,
2290 allow_cf, cq_s, perf_press, perf_rates, deferred_logger);
2291 for (int comp = 0; comp < this->num_conservation_quantities_; ++comp) {
2292 well_q_s[comp] += cq_s[comp];
2293 }
2294 }
2295 }
2296 const auto& comm = this->parallel_well_info_.communication();
2297 if (comm.size() > 1)
2298 {
2299 comm.sum(well_q_s.data(), well_q_s.size());
2300 }
2301 return well_q_s;
2302 }
2303
2304
2305 template <typename TypeTag>
2306 std::vector<typename MultisegmentWell<TypeTag>::Scalar>
2308 getPrimaryVars() const
2309 {
2310 const int num_seg = this->numberOfSegments();
2311 constexpr int num_eq = MSWEval::numWellEq;
2312 std::vector<Scalar> retval(num_seg * num_eq);
2313 for (int ii = 0; ii < num_seg; ++ii) {
2314 const auto& pv = this->primary_variables_.value(ii);
2315 std::copy(pv.begin(), pv.end(), retval.begin() + ii * num_eq);
2316 }
2317 return retval;
2318 }
2319
2320
2321
2322
2323 template <typename TypeTag>
2324 int
2326 setPrimaryVars(typename std::vector<Scalar>::const_iterator it)
2327 {
2328 const int num_seg = this->numberOfSegments();
2329 constexpr int num_eq = MSWEval::numWellEq;
2330 std::array<Scalar, num_eq> tmp;
2331 for (int ii = 0; ii < num_seg; ++ii) {
2332 const auto start = it + ii * num_eq;
2333 std::copy(start, start + num_eq, tmp.begin());
2334 this->primary_variables_.setValue(ii, tmp);
2335 }
2336 return num_seg * num_eq;
2337 }
2338
2339 template <typename TypeTag>
2342 getFirstPerforationFluidStateInfo(const Simulator& simulator) const
2343 {
2344 Scalar fsTemperature = 0.0;
2345 using SaltConcType = typename std::decay<decltype(std::declval<decltype(simulator.model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type;
2346 SaltConcType fsSaltConcentration{};
2347
2348 // If this process does not contain active perforations, this->well_cells_ is empty.
2349 if (this->well_cells_.size() > 0) {
2350 // We use the pvt region of first perforated cell, so we look for global index 0
2351 // TODO: it should be a member of the WellInterface, initialized properly
2352 const int cell_idx = this->well_cells_[0];
2353 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
2354 const auto& fs = intQuants.fluidState();
2355
2356 fsTemperature = fs.temperature(FluidSystem::oilPhaseIdx).value();
2357 fsSaltConcentration = fs.saltConcentration();
2358 }
2359
2360 auto info = std::make_tuple(fsTemperature, fsSaltConcentration);
2361
2362 // The following broadcast call is neccessary to ensure that processes that do *not* contain
2363 // the first perforation get the correct temperature, saltConcentration and pvt_region_index
2364 return this->parallel_well_info_.communication().size() == 1 ? info : this->pw_info_.broadcastFirstPerforationValue(info);
2365 }
2366
2367} // namespace Opm
2368
2369#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:2102
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellStateType &well_state, DeferredLogger &deferred_logger, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:719
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: MultisegmentWell_impl.hpp:2069
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:891
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:856
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:878
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:1521
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_indx) const
Definition: MultisegmentWell_impl.hpp:1164
void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:614
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1224
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:2122
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2308
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=false, const bool fixed_status=false) override
Definition: MultisegmentWell_impl.hpp:1646
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
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:2262
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:1232
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1184
void updateIPRImplicit(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1399
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2007
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1133
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:2326
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:2016
int debug_cost_counter_
Definition: MultisegmentWell.hpp:182
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:774
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:1060
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:2179
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:643
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:1298
FSInfo getFirstPerforationFluidStateInfo(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2342
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
void assembleWellEqWithoutIteration(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, const GroupState< Scalar > &group_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1824
Scalar maxPerfPress(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2234
void computeInitialSegmentFluids(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:700
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1468
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:2080
void calculateExplicitQuantities(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:758
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 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:2049
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:2062
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:116
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:580
bool thp_update_iterations
Definition: WellInterface.hpp:376
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