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