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