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