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