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