StandardWell_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4 Copyright 2016 - 2017 IRIS AS.
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20*/
21
22// Improve IDE experience
23#ifndef OPM_STANDARDWELL_HEADER_INCLUDED
24#define OPM_STANDARDWELL_IMPL_HEADER_INCLUDED
25#include <config.h>
27#endif
28
29#include <opm/common/Exceptions.hpp>
30
31#include <opm/input/eclipse/Units/Units.hpp>
32
33#include <opm/material/densead/EvaluationFormat.hpp>
34
40
41#include <fmt/format.h>
42
43#include <algorithm>
44#include <cstddef>
45#include <functional>
46#include <numeric>
47
48namespace {
49
50template<class dValue, class Value>
51auto dValueError(const dValue& d,
52 const std::string& name,
53 const std::string& methodName,
54 const Value& Rs,
55 const Value& Rv,
56 const Value& pressure)
57{
58 return fmt::format("Problematic d value {} obtained for well {}"
59 " during {} calculations with rs {}"
60 ", rv {} and pressure {}."
61 " Continue as if no dissolution (rs = 0) and vaporization (rv = 0) "
62 " for this connection.", d, name, methodName, Rs, Rv, pressure);
63}
64
65}
66
67namespace Opm
68{
69
70 template<typename TypeTag>
72 StandardWell(const Well& well,
73 const ParallelWellInfo<Scalar>& pw_info,
74 const int time_step,
75 const ModelParameters& param,
76 const RateConverterType& rate_converter,
77 const int pvtRegionIdx,
78 const int num_components,
79 const int num_phases,
80 const int index_of_well,
81 const std::vector<PerforationData<Scalar>>& perf_data)
82 : Base(well, pw_info, time_step, param, rate_converter, pvtRegionIdx, num_components, num_phases, index_of_well, perf_data)
83 , StdWellEval(static_cast<const WellInterfaceIndices<FluidSystem,Indices>&>(*this))
84 , regularize_(false)
85 {
87 }
88
89
90
91
92
93 template<typename TypeTag>
94 void
96 init(const PhaseUsage* phase_usage_arg,
97 const std::vector<Scalar>& depth_arg,
98 const Scalar gravity_arg,
99 const int num_cells,
100 const std::vector< Scalar >& B_avg,
101 const bool changed_to_open_this_step)
102 {
103 Base::init(phase_usage_arg, depth_arg, gravity_arg, num_cells, B_avg, changed_to_open_this_step);
104 this->StdWellEval::init(this->perf_depth_, depth_arg, num_cells, Base::has_polymermw);
105 }
106
107
108
109
110
111 template<typename TypeTag>
114 {
115 this->primary_variables_.init();
116 }
117
118
119
120
121
122 template<typename TypeTag>
123 template<class Value>
124 void
127 const std::vector<Value>& mob,
128 const Value& bhp,
129 const std::vector<Scalar>& Tw,
130 const int perf,
131 const bool allow_cf,
132 std::vector<Value>& cq_s,
133 PerforationRates<Scalar>& perf_rates,
134 DeferredLogger& deferred_logger) const
135 {
136 auto obtain = [this](const Eval& value)
137 {
138 if constexpr (std::is_same_v<Value, Scalar>) {
139 static_cast<void>(this); // suppress clang warning
140 return getValue(value);
141 } else {
142 return this->extendEval(value);
143 }
144 };
145 auto obtainN = [](const auto& value)
146 {
147 if constexpr (std::is_same_v<Value, Scalar>) {
148 return getValue(value);
149 } else {
150 return value;
151 }
152 };
153 auto zeroElem = [this]()
154 {
155 if constexpr (std::is_same_v<Value, Scalar>) {
156 static_cast<void>(this); // suppress clang warning
157 return 0.0;
158 } else {
159 return Value{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
160 }
161 };
162
163 const auto& fs = intQuants.fluidState();
164 const Value pressure = obtain(this->getPerfCellPressure(fs));
165 const Value rs = obtain(fs.Rs());
166 const Value rv = obtain(fs.Rv());
167 const Value rvw = obtain(fs.Rvw());
168 const Value rsw = obtain(fs.Rsw());
169
170 std::vector<Value> b_perfcells_dense(this->numComponents(), zeroElem());
171 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
172 if (!FluidSystem::phaseIsActive(phaseIdx)) {
173 continue;
174 }
175 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
176 b_perfcells_dense[compIdx] = obtain(fs.invB(phaseIdx));
177 }
178 if constexpr (has_solvent) {
179 b_perfcells_dense[Indices::contiSolventEqIdx] = obtain(intQuants.solventInverseFormationVolumeFactor());
180 }
181
182 if constexpr (has_zFraction) {
183 if (this->isInjector()) {
184 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
185 b_perfcells_dense[gasCompIdx] *= (1.0 - this->wsolvent());
186 b_perfcells_dense[gasCompIdx] += this->wsolvent()*intQuants.zPureInvFormationVolumeFactor().value();
187 }
188 }
189
190 Value skin_pressure = zeroElem();
191 if (has_polymermw) {
192 if (this->isInjector()) {
193 const int pskin_index = Bhp + 1 + this->numPerfs() + perf;
194 skin_pressure = obtainN(this->primary_variables_.eval(pskin_index));
195 }
196 }
197
198 // surface volume fraction of fluids within wellbore
199 std::vector<Value> cmix_s(this->numComponents(), zeroElem());
200 for (int componentIdx = 0; componentIdx < this->numComponents(); ++componentIdx) {
201 cmix_s[componentIdx] = obtainN(this->primary_variables_.surfaceVolumeFraction(componentIdx));
202 }
203
204 computePerfRate(mob,
205 pressure,
206 bhp,
207 rs,
208 rv,
209 rvw,
210 rsw,
211 b_perfcells_dense,
212 Tw,
213 perf,
214 allow_cf,
215 skin_pressure,
216 cmix_s,
217 cq_s,
218 perf_rates,
219 deferred_logger);
220 }
221
222 template<typename TypeTag>
223 template<class Value>
224 void
226 computePerfRate(const std::vector<Value>& mob,
227 const Value& pressure,
228 const Value& bhp,
229 const Value& rs,
230 const Value& rv,
231 const Value& rvw,
232 const Value& rsw,
233 std::vector<Value>& b_perfcells_dense,
234 const std::vector<Scalar>& Tw,
235 const int perf,
236 const bool allow_cf,
237 const Value& skin_pressure,
238 const std::vector<Value>& cmix_s,
239 std::vector<Value>& cq_s,
240 PerforationRates<Scalar>& perf_rates,
241 DeferredLogger& deferred_logger) const
242 {
243 // Pressure drawdown (also used to determine direction of flow)
244 const Value well_pressure = bhp + this->connections_.pressure_diff(perf);
245 Value drawdown = pressure - well_pressure;
246 if (this->isInjector()) {
247 drawdown += skin_pressure;
248 }
249
250 // producing perforations
251 if (drawdown > 0) {
252 // Do nothing if crossflow is not allowed
253 if (!allow_cf && this->isInjector()) {
254 return;
255 }
256
257 // compute component volumetric rates at standard conditions
258 for (int componentIdx = 0; componentIdx < this->numComponents(); ++componentIdx) {
259 const Value cq_p = - Tw[componentIdx] * (mob[componentIdx] * drawdown);
260 cq_s[componentIdx] = b_perfcells_dense[componentIdx] * cq_p;
261 }
262
263 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
264 gasOilPerfRateProd(cq_s, perf_rates, rv, rs, rvw);
265 } else if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
266 gasWaterPerfRateProd(cq_s, perf_rates, rvw, rsw);
267 }
268 } else {
269 // Do nothing if crossflow is not allowed
270 if (!allow_cf && this->isProducer()) {
271 return;
272 }
273
274 // Using total mobilities
275 Value total_mob_dense = mob[0];
276 for (int componentIdx = 1; componentIdx < this->numComponents(); ++componentIdx) {
277 total_mob_dense += mob[componentIdx];
278 }
279
280 // compute volume ratio between connection at standard conditions
281 Value volumeRatio = bhp * 0.0; // initialize it with the correct type
282
283 if (FluidSystem::enableVaporizedWater() && FluidSystem::enableDissolvedGasInWater()) {
284 disOilVapWatVolumeRatio(volumeRatio, rvw, rsw, pressure,
285 cmix_s, b_perfcells_dense, deferred_logger);
286 // DISGASW only supported for gas-water CO2STORE/H2STORE case
287 // and the simulator will throw long before it reach to this point in the code
288 // For blackoil support of DISGASW we need to add the oil component here
289 assert(FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx));
290 assert(FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx));
291 assert(!FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx));
292 } else {
293
294 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
295 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
296 volumeRatio += cmix_s[waterCompIdx] / b_perfcells_dense[waterCompIdx];
297 }
298
299 if constexpr (Indices::enableSolvent) {
300 volumeRatio += cmix_s[Indices::contiSolventEqIdx] / b_perfcells_dense[Indices::contiSolventEqIdx];
301 }
302
303 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
304 gasOilVolumeRatio(volumeRatio, rv, rs, pressure,
305 cmix_s, b_perfcells_dense, deferred_logger);
306 } else {
307 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
308 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
309 volumeRatio += cmix_s[oilCompIdx] / b_perfcells_dense[oilCompIdx];
310 }
311 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
312 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
313 volumeRatio += cmix_s[gasCompIdx] / b_perfcells_dense[gasCompIdx];
314 }
315 }
316 }
317
318 // injecting connections total volumerates at standard conditions
319 for (int componentIdx = 0; componentIdx < this->numComponents(); ++componentIdx) {
320 const Value cqt_i = - Tw[componentIdx] * (total_mob_dense * drawdown);
321 Value cqt_is = cqt_i / volumeRatio;
322 cq_s[componentIdx] = cmix_s[componentIdx] * cqt_is;
323 }
324
325 // calculating the perforation solution gas rate and solution oil rates
326 if (this->isProducer()) {
327 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
328 gasOilPerfRateInj(cq_s, perf_rates,
329 rv, rs, pressure, rvw, deferred_logger);
330 }
331 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
332 //no oil
333 gasWaterPerfRateInj(cq_s, perf_rates, rvw, rsw,
334 pressure, deferred_logger);
335 }
336 }
337 }
338 }
339
340
341 template<typename TypeTag>
342 void
345 const double dt,
346 const Well::InjectionControls& inj_controls,
347 const Well::ProductionControls& prod_controls,
348 WellState<Scalar>& well_state,
349 const GroupState<Scalar>& group_state,
350 DeferredLogger& deferred_logger)
351 {
352 // TODO: only_wells should be put back to save some computation
353 // for example, the matrices B C does not need to update if only_wells
354 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
355
356 // clear all entries
357 this->linSys_.clear();
358
359 assembleWellEqWithoutIterationImpl(simulator, dt, inj_controls,
360 prod_controls, well_state,
361 group_state, deferred_logger);
362 }
363
364
365
366
367 template<typename TypeTag>
368 void
371 const double dt,
372 const Well::InjectionControls& inj_controls,
373 const Well::ProductionControls& prod_controls,
374 WellState<Scalar>& well_state,
375 const GroupState<Scalar>& group_state,
376 DeferredLogger& deferred_logger)
377 {
378 // try to regularize equation if the well does not converge
379 const Scalar regularization_factor = this->regularize_? this->param_.regularization_factor_wells_ : 1.0;
380 const Scalar volume = 0.1 * unit::cubic(unit::feet) * regularization_factor;
381
382 auto& ws = well_state.well(this->index_of_well_);
383 ws.phase_mixing_rates.fill(0.0);
384
385
386 const int np = this->number_of_phases_;
387
388 std::vector<RateVector> connectionRates = this->connectionRates_; // Copy to get right size.
389
390 auto& perf_data = ws.perf_data;
391 auto& perf_rates = perf_data.phase_rates;
392 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
393 // Calculate perforation quantities.
394 std::vector<EvalWell> cq_s(this->num_components_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.0});
395 EvalWell water_flux_s{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
396 EvalWell cq_s_zfrac_effective{this->primary_variables_.numWellEq() + Indices::numEq, 0.0};
397 calculateSinglePerf(simulator, perf, well_state, connectionRates,
398 cq_s, water_flux_s, cq_s_zfrac_effective, deferred_logger);
399
400 // Equation assembly for this perforation.
401 if constexpr (has_polymer && Base::has_polymermw) {
402 if (this->isInjector()) {
403 handleInjectivityEquations(simulator, well_state, perf,
404 water_flux_s, deferred_logger);
405 }
406 }
407 const int cell_idx = this->well_cells_[perf];
408 for (int componentIdx = 0; componentIdx < this->num_components_; ++componentIdx) {
409 // the cq_s entering mass balance equations need to consider the efficiency factors.
410 const EvalWell cq_s_effective = cq_s[componentIdx] * this->well_efficiency_factor_;
411
412 connectionRates[perf][componentIdx] = Base::restrictEval(cq_s_effective);
413
415 assemblePerforationEq(cq_s_effective,
416 componentIdx,
417 cell_idx,
418 this->primary_variables_.numWellEq(),
419 this->linSys_);
420
421 // Store the perforation phase flux for later usage.
422 if (has_solvent && componentIdx == Indices::contiSolventEqIdx) {
423 auto& perf_rate_solvent = perf_data.solvent_rates;
424 perf_rate_solvent[perf] = cq_s[componentIdx].value();
425 } else {
426 perf_rates[perf*np + this->modelCompIdxToFlowCompIdx(componentIdx)] = cq_s[componentIdx].value();
427 }
428 }
429
430 if constexpr (has_zFraction) {
432 assembleZFracEq(cq_s_zfrac_effective,
433 cell_idx,
434 this->primary_variables_.numWellEq(),
435 this->linSys_);
436 }
437 }
438 // Update the connection
439 this->connectionRates_ = connectionRates;
440
441 // Accumulate dissolved gas and vaporized oil flow rates across all
442 // ranks sharing this well (this->index_of_well_).
443 {
444 const auto& comm = this->parallel_well_info_.communication();
445 comm.sum(ws.phase_mixing_rates.data(), ws.phase_mixing_rates.size());
446 }
447
448 // accumulate resWell_ and duneD_ in parallel to get effects of all perforations (might be distributed)
449 this->linSys_.sumDistributed(this->parallel_well_info_.communication());
450
451 // add vol * dF/dt + Q to the well equations;
452 for (int componentIdx = 0; componentIdx < numWellConservationEq; ++componentIdx) {
453 // TODO: following the development in MSW, we need to convert the volume of the wellbore to be surface volume
454 // since all the rates are under surface condition
455 EvalWell resWell_loc(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
456 if (FluidSystem::numActivePhases() > 1) {
457 assert(dt > 0);
458 resWell_loc += (this->primary_variables_.surfaceVolumeFraction(componentIdx) -
459 this->F0_[componentIdx]) * volume / dt;
460 }
461 resWell_loc -= this->primary_variables_.getQs(componentIdx) * this->well_efficiency_factor_;
463 assembleSourceEq(resWell_loc,
464 componentIdx,
465 this->primary_variables_.numWellEq(),
466 this->linSys_);
467 }
468
469 const auto& summaryState = simulator.vanguard().summaryState();
470 const Schedule& schedule = simulator.vanguard().schedule();
471 const bool stopped_or_zero_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
473 assembleControlEq(well_state, group_state,
474 schedule, summaryState,
475 inj_controls, prod_controls,
476 this->primary_variables_,
477 this->connections_.rho(),
478 this->linSys_,
479 stopped_or_zero_target,
480 deferred_logger);
481
482
483 // do the local inversion of D.
484 try {
485 this->linSys_.invert();
486 } catch( ... ) {
487 OPM_DEFLOG_PROBLEM(NumericalProblem, "Error when inverting local well equations for well " + name(), deferred_logger);
488 }
489 }
490
491
492
493
494 template<typename TypeTag>
495 void
497 calculateSinglePerf(const Simulator& simulator,
498 const int perf,
499 WellState<Scalar>& well_state,
500 std::vector<RateVector>& connectionRates,
501 std::vector<EvalWell>& cq_s,
502 EvalWell& water_flux_s,
503 EvalWell& cq_s_zfrac_effective,
504 DeferredLogger& deferred_logger) const
505 {
506 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
507 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
508 const int cell_idx = this->well_cells_[perf];
509 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
510 std::vector<EvalWell> mob(this->num_components_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
511 getMobility(simulator, perf, mob, deferred_logger);
512
513 PerforationRates<Scalar> perf_rates;
514 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
515 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
516 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
517 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
518 cq_s, perf_rates, deferred_logger);
519
520 auto& ws = well_state.well(this->index_of_well_);
521 auto& perf_data = ws.perf_data;
522 if constexpr (has_polymer && Base::has_polymermw) {
523 if (this->isInjector()) {
524 // Store the original water flux computed from the reservoir quantities.
525 // It will be required to assemble the injectivity equations.
526 const unsigned water_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
527 water_flux_s = cq_s[water_comp_idx];
528 // Modify the water flux for the rest of this function to depend directly on the
529 // local water velocity primary variable.
530 handleInjectivityRate(simulator, perf, cq_s);
531 }
532 }
533
534 // updating the solution gas rate and solution oil rate
535 if (this->isProducer()) {
536 ws.phase_mixing_rates[ws.dissolved_gas] += perf_rates.dis_gas;
537 ws.phase_mixing_rates[ws.dissolved_gas_in_water] += perf_rates.dis_gas_in_water;
538 ws.phase_mixing_rates[ws.vaporized_oil] += perf_rates.vap_oil;
539 ws.phase_mixing_rates[ws.vaporized_water] += perf_rates.vap_wat;
540 perf_data.phase_mixing_rates[perf][ws.dissolved_gas] = perf_rates.dis_gas;
541 perf_data.phase_mixing_rates[perf][ws.dissolved_gas_in_water] = perf_rates.dis_gas_in_water;
542 perf_data.phase_mixing_rates[perf][ws.vaporized_oil] = perf_rates.vap_oil;
543 perf_data.phase_mixing_rates[perf][ws.vaporized_water] = perf_rates.vap_wat;
544 }
545
546 if constexpr (has_energy) {
547 connectionRates[perf][Indices::contiEnergyEqIdx] =
548 connectionRateEnergy(simulator.problem().maxOilSaturation(cell_idx),
549 cq_s, intQuants, deferred_logger);
550 }
551
552 if constexpr (has_polymer) {
553 std::variant<Scalar,EvalWell> polymerConcentration;
554 if (this->isInjector()) {
555 polymerConcentration = this->wpolymer();
556 } else {
557 polymerConcentration = this->extendEval(intQuants.polymerConcentration() *
558 intQuants.polymerViscosityCorrection());
559 }
560
561 [[maybe_unused]] EvalWell cq_s_poly;
562 std::tie(connectionRates[perf][Indices::contiPolymerEqIdx],
563 cq_s_poly) =
564 this->connections_.connectionRatePolymer(perf_data.polymer_rates[perf],
565 cq_s, polymerConcentration);
566
567 if constexpr (Base::has_polymermw) {
568 updateConnectionRatePolyMW(cq_s_poly, intQuants, well_state,
569 perf, connectionRates, deferred_logger);
570 }
571 }
572
573 if constexpr (has_foam) {
574 std::variant<Scalar,EvalWell> foamConcentration;
575 if (this->isInjector()) {
576 foamConcentration = this->wfoam();
577 } else {
578 foamConcentration = this->extendEval(intQuants.foamConcentration());
579 }
580 connectionRates[perf][Indices::contiFoamEqIdx] =
581 this->connections_.connectionRateFoam(cq_s, foamConcentration,
582 FoamModule::transportPhase(),
583 deferred_logger);
584 }
585
586 if constexpr (has_zFraction) {
587 std::variant<Scalar,std::array<EvalWell,2>> solventConcentration;
588 if (this->isInjector()) {
589 solventConcentration = this->wsolvent();
590 } else {
591 solventConcentration = std::array{this->extendEval(intQuants.xVolume()),
592 this->extendEval(intQuants.yVolume())};
593 }
594 std::tie(connectionRates[perf][Indices::contiZfracEqIdx],
595 cq_s_zfrac_effective) =
596 this->connections_.connectionRatezFraction(perf_data.solvent_rates[perf],
597 perf_rates.dis_gas, cq_s,
598 solventConcentration);
599 }
600
601 if constexpr (has_brine) {
602 std::variant<Scalar,EvalWell> saltConcentration;
603 if (this->isInjector()) {
604 saltConcentration = this->wsalt();
605 } else {
606 saltConcentration = this->extendEval(intQuants.fluidState().saltConcentration());
607 }
608
609 connectionRates[perf][Indices::contiBrineEqIdx] =
610 this->connections_.connectionRateBrine(perf_data.brine_rates[perf],
611 perf_rates.vap_wat, cq_s,
612 saltConcentration);
613 }
614
615 if constexpr (has_micp) {
616 std::variant<Scalar,EvalWell> microbialConcentration;
617 std::variant<Scalar,EvalWell> oxygenConcentration;
618 std::variant<Scalar,EvalWell> ureaConcentration;
619 if (this->isInjector()) {
620 microbialConcentration = this->wmicrobes();
621 oxygenConcentration = this->woxygen();
622 ureaConcentration = this->wurea();
623 } else {
624 microbialConcentration = this->extendEval(intQuants.microbialConcentration());
625 oxygenConcentration = this->extendEval(intQuants.oxygenConcentration());
626 ureaConcentration = this->extendEval(intQuants.ureaConcentration());
627 }
628 std::tie(connectionRates[perf][Indices::contiMicrobialEqIdx],
629 connectionRates[perf][Indices::contiOxygenEqIdx],
630 connectionRates[perf][Indices::contiUreaEqIdx]) =
631 this->connections_.connectionRatesMICP(cq_s,
632 microbialConcentration,
633 oxygenConcentration,
634 ureaConcentration);
635 }
636
637 // Store the perforation pressure for later usage.
638 perf_data.pressure[perf] = ws.bhp + this->connections_.pressure_diff(perf);
639 }
640
641
642
643 template<typename TypeTag>
644 template<class Value>
645 void
647 getMobility(const Simulator& simulator,
648 const int perf,
649 std::vector<Value>& mob,
650 DeferredLogger& deferred_logger) const
651 {
652 auto obtain = [this](const Eval& value)
653 {
654 if constexpr (std::is_same_v<Value, Scalar>) {
655 static_cast<void>(this); // suppress clang warning
656 return getValue(value);
657 } else {
658 return this->extendEval(value);
659 }
660 };
661 WellInterface<TypeTag>::getMobility(simulator, perf, mob,
662 obtain, deferred_logger);
663
664 // modify the water mobility if polymer is present
665 if constexpr (has_polymer) {
666 if (!FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
667 OPM_DEFLOG_THROW(std::runtime_error, "Water is required when polymer is active", deferred_logger);
668 }
669
670 // for the cases related to polymer molecular weight, we assume fully mixing
671 // as a result, the polymer and water share the same viscosity
672 if constexpr (!Base::has_polymermw) {
673 if constexpr (std::is_same_v<Value, Scalar>) {
674 std::vector<EvalWell> mob_eval(this->num_components_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
675 for (std::size_t i = 0; i < mob.size(); ++i) {
676 mob_eval[i].setValue(mob[i]);
677 }
678 updateWaterMobilityWithPolymer(simulator, perf, mob_eval, deferred_logger);
679 for (std::size_t i = 0; i < mob.size(); ++i) {
680 mob[i] = getValue(mob_eval[i]);
681 }
682 } else {
683 updateWaterMobilityWithPolymer(simulator, perf, mob, deferred_logger);
684 }
685 }
686 }
687
688 // if the injecting well has WINJMULT setup, we update the mobility accordingly
689 if (this->isInjector() && this->well_ecl_.getInjMultMode() != Well::InjMultMode::NONE) {
690 const Scalar bhp = this->primary_variables_.value(Bhp);
691 const Scalar perf_press = bhp + this->connections_.pressure_diff(perf);
692 const Scalar multiplier = this->getInjMult(perf, bhp, perf_press);
693 for (std::size_t i = 0; i < mob.size(); ++i) {
694 mob[i] *= multiplier;
695 }
696 }
697 }
698
699
700 template<typename TypeTag>
701 void
703 updateWellState(const Simulator& simulator,
704 const BVectorWell& dwells,
705 WellState<Scalar>& well_state,
706 DeferredLogger& deferred_logger)
707 {
708 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
709
710 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
711 updatePrimaryVariablesNewton(dwells, stop_or_zero_rate_target, deferred_logger);
712
713 const auto& summary_state = simulator.vanguard().summaryState();
714 updateWellStateFromPrimaryVariables(stop_or_zero_rate_target, well_state, summary_state, deferred_logger);
715 Base::calculateReservoirRates(simulator.vanguard().eclState().runspec().co2Storage(), well_state.well(this->index_of_well_));
716 }
717
718
719
720
721
722 template<typename TypeTag>
723 void
726 const bool stop_or_zero_rate_target,
727 DeferredLogger& deferred_logger)
728 {
729 const Scalar dFLimit = this->param_.dwell_fraction_max_;
730 const Scalar dBHPLimit = this->param_.dbhp_max_rel_;
731 this->primary_variables_.updateNewton(dwells, stop_or_zero_rate_target, dFLimit, dBHPLimit, deferred_logger);
732
733 // for the water velocity and skin pressure
734 if constexpr (Base::has_polymermw) {
735 this->primary_variables_.updateNewtonPolyMW(dwells);
736 }
737
738 this->primary_variables_.checkFinite(deferred_logger);
739 }
740
741
742
743
744
745 template<typename TypeTag>
746 void
748 updateWellStateFromPrimaryVariables(const bool stop_or_zero_rate_target,
749 WellState<Scalar>& well_state,
750 const SummaryState& summary_state,
751 DeferredLogger& deferred_logger) const
752 {
753 this->StdWellEval::updateWellStateFromPrimaryVariables(stop_or_zero_rate_target, well_state, summary_state, deferred_logger);
754
755 // other primary variables related to polymer injectivity study
756 if constexpr (Base::has_polymermw) {
757 this->primary_variables_.copyToWellStatePolyMW(well_state);
758 }
759 }
760
761
762
763
764
765 template<typename TypeTag>
766 void
768 updateIPR(const Simulator& simulator, DeferredLogger& deferred_logger) const
769 {
770 // TODO: not handling solvent related here for now
771
772 // initialize all the values to be zero to begin with
773 std::fill(this->ipr_a_.begin(), this->ipr_a_.end(), 0.);
774 std::fill(this->ipr_b_.begin(), this->ipr_b_.end(), 0.);
775
776 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
777 std::vector<Scalar> mob(this->num_components_, 0.0);
778 getMobility(simulator, perf, mob, deferred_logger);
779
780 const int cell_idx = this->well_cells_[perf];
781 const auto& int_quantities = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
782 const auto& fs = int_quantities.fluidState();
783 // the pressure of the reservoir grid block the well connection is in
784 Scalar p_r = this->getPerfCellPressure(fs).value();
785
786 // calculating the b for the connection
787 std::vector<Scalar> b_perf(this->num_components_);
788 for (std::size_t phase = 0; phase < FluidSystem::numPhases; ++phase) {
789 if (!FluidSystem::phaseIsActive(phase)) {
790 continue;
791 }
792 const unsigned comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phase));
793 b_perf[comp_idx] = fs.invB(phase).value();
794 }
795 if constexpr (has_solvent) {
796 b_perf[Indices::contiSolventEqIdx] = int_quantities.solventInverseFormationVolumeFactor().value();
797 }
798
799 // the pressure difference between the connection and BHP
800 const Scalar h_perf = this->connections_.pressure_diff(perf);
801 const Scalar pressure_diff = p_r - h_perf;
802
803 // Let us add a check, since the pressure is calculated based on zero value BHP
804 // it should not be negative anyway. If it is negative, we might need to re-formulate
805 // to taking into consideration the crossflow here.
806 if ( (this->isProducer() && pressure_diff < 0.) || (this->isInjector() && pressure_diff > 0.) ) {
807 deferred_logger.debug("CROSSFLOW_IPR",
808 "cross flow found when updateIPR for well " + name()
809 + " . The connection is ignored in IPR calculations");
810 // we ignore these connections for now
811 continue;
812 }
813
814 // the well index associated with the connection
815 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quantities, cell_idx);
816 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
817 const std::vector<Scalar> tw_perf = this->wellIndex(perf,
818 int_quantities,
819 trans_mult,
820 wellstate_nupcol);
821 std::vector<Scalar> ipr_a_perf(this->ipr_a_.size());
822 std::vector<Scalar> ipr_b_perf(this->ipr_b_.size());
823 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
824 const Scalar tw_mob = tw_perf[comp_idx] * mob[comp_idx] * b_perf[comp_idx];
825 ipr_a_perf[comp_idx] += tw_mob * pressure_diff;
826 ipr_b_perf[comp_idx] += tw_mob;
827 }
828
829 // we need to handle the rs and rv when both oil and gas are present
830 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
831 const unsigned oil_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
832 const unsigned gas_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
833 const Scalar rs = (fs.Rs()).value();
834 const Scalar rv = (fs.Rv()).value();
835
836 const Scalar dis_gas_a = rs * ipr_a_perf[oil_comp_idx];
837 const Scalar vap_oil_a = rv * ipr_a_perf[gas_comp_idx];
838
839 ipr_a_perf[gas_comp_idx] += dis_gas_a;
840 ipr_a_perf[oil_comp_idx] += vap_oil_a;
841
842 const Scalar dis_gas_b = rs * ipr_b_perf[oil_comp_idx];
843 const Scalar vap_oil_b = rv * ipr_b_perf[gas_comp_idx];
844
845 ipr_b_perf[gas_comp_idx] += dis_gas_b;
846 ipr_b_perf[oil_comp_idx] += vap_oil_b;
847 }
848
849 for (std::size_t comp_idx = 0; comp_idx < ipr_a_perf.size(); ++comp_idx) {
850 this->ipr_a_[comp_idx] += ipr_a_perf[comp_idx];
851 this->ipr_b_[comp_idx] += ipr_b_perf[comp_idx];
852 }
853 }
854 this->parallel_well_info_.communication().sum(this->ipr_a_.data(), this->ipr_a_.size());
855 this->parallel_well_info_.communication().sum(this->ipr_b_.data(), this->ipr_b_.size());
856 }
857
858 template<typename TypeTag>
859 void
861 updateIPRImplicit(const Simulator& simulator,
862 WellState<Scalar>& well_state,
863 DeferredLogger& deferred_logger)
864 {
865 // Compute IPR based on *converged* well-equation:
866 // For a component rate r the derivative dr/dbhp is obtained by
867 // dr/dbhp = - (partial r/partial x) * inv(partial Eq/partial x) * (partial Eq/partial bhp_target)
868 // where Eq(x)=0 is the well equation setup with bhp control and primary variables x
869
870 // We shouldn't have zero rates at this stage, but check
871 bool zero_rates;
872 auto rates = well_state.well(this->index_of_well_).surface_rates;
873 zero_rates = true;
874 for (std::size_t p = 0; p < rates.size(); ++p) {
875 zero_rates &= rates[p] == 0.0;
876 }
877 auto& ws = well_state.well(this->index_of_well_);
878 if (zero_rates) {
879 const auto msg = fmt::format("updateIPRImplicit: Well {} has zero rate, IPRs might be problematic", this->name());
880 deferred_logger.debug(msg);
881 /*
882 // could revert to standard approach here:
883 updateIPR(simulator, deferred_logger);
884 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx){
885 const int idx = this->modelCompIdxToFlowCompIdx(comp_idx);
886 ws.implicit_ipr_a[idx] = this->ipr_a_[comp_idx];
887 ws.implicit_ipr_b[idx] = this->ipr_b_[comp_idx];
888 }
889 return;
890 */
891 }
892 const auto& group_state = simulator.problem().wellModel().groupState();
893
894 std::fill(ws.implicit_ipr_a.begin(), ws.implicit_ipr_a.end(), 0.);
895 std::fill(ws.implicit_ipr_b.begin(), ws.implicit_ipr_b.end(), 0.);
896
897 auto inj_controls = Well::InjectionControls(0);
898 auto prod_controls = Well::ProductionControls(0);
899 prod_controls.addControl(Well::ProducerCMode::BHP);
900 prod_controls.bhp_limit = well_state.well(this->index_of_well_).bhp;
901
902 // Set current control to bhp, and bhp value in state, modify bhp limit in control object.
903 const auto cmode = ws.production_cmode;
904 ws.production_cmode = Well::ProducerCMode::BHP;
905 const double dt = simulator.timeStepSize();
906 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
907
908 const size_t nEq = this->primary_variables_.numWellEq();
909 BVectorWell rhs(1);
910 rhs[0].resize(nEq);
911 // rhs = 0 except -1 for control eq
912 for (size_t i=0; i < nEq; ++i){
913 rhs[0][i] = 0.0;
914 }
915 rhs[0][Bhp] = -1.0;
916
917 BVectorWell x_well(1);
918 x_well[0].resize(nEq);
919 this->linSys_.solve(rhs, x_well);
920
921 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx){
922 EvalWell comp_rate = this->primary_variables_.getQs(comp_idx);
923 const int idx = this->modelCompIdxToFlowCompIdx(comp_idx);
924 for (size_t pvIdx = 0; pvIdx < nEq; ++pvIdx) {
925 // well primary variable derivatives in EvalWell start at position Indices::numEq
926 ws.implicit_ipr_b[idx] -= x_well[0][pvIdx]*comp_rate.derivative(pvIdx+Indices::numEq);
927 }
928 ws.implicit_ipr_a[idx] = ws.implicit_ipr_b[idx]*ws.bhp - comp_rate.value();
929 }
930 // reset cmode
931 ws.production_cmode = cmode;
932 }
933
934 template<typename TypeTag>
935 void
938 const Simulator& simulator,
939 DeferredLogger& deferred_logger)
940 {
941 const auto& summaryState = simulator.vanguard().summaryState();
942 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
943 // Crude but works: default is one atmosphere.
944 // TODO: a better way to detect whether the BHP is defaulted or not
945 const bool bhp_limit_not_defaulted = bhp_limit > 1.5 * unit::barsa;
946 if ( bhp_limit_not_defaulted || !this->wellHasTHPConstraints(summaryState) ) {
947 // if the BHP limit is not defaulted or the well does not have a THP limit
948 // we need to check the BHP limit
949 Scalar total_ipr_mass_rate = 0.0;
950 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx)
951 {
952 if (!FluidSystem::phaseIsActive(phaseIdx)) {
953 continue;
954 }
955
956 const unsigned compIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
957 const Scalar ipr_rate = this->ipr_a_[compIdx] - this->ipr_b_[compIdx] * bhp_limit;
958
959 const Scalar rho = FluidSystem::referenceDensity( phaseIdx, Base::pvtRegionIdx() );
960 total_ipr_mass_rate += ipr_rate * rho;
961 }
962 if ( (this->isProducer() && total_ipr_mass_rate < 0.) || (this->isInjector() && total_ipr_mass_rate > 0.) ) {
963 this->operability_status_.operable_under_only_bhp_limit = false;
964 }
965
966 // checking whether running under BHP limit will violate THP limit
967 if (this->operability_status_.operable_under_only_bhp_limit && this->wellHasTHPConstraints(summaryState)) {
968 // option 1: calculate well rates based on the BHP limit.
969 // option 2: stick with the above IPR curve
970 // we use IPR here
971 std::vector<Scalar> well_rates_bhp_limit;
972 computeWellRatesWithBhp(simulator, bhp_limit, well_rates_bhp_limit, deferred_logger);
973
974 this->adaptRatesForVFP(well_rates_bhp_limit);
975 const Scalar thp_limit = this->getTHPConstraint(summaryState);
976 const Scalar thp = WellBhpThpCalculator(*this).calculateThpFromBhp(well_rates_bhp_limit,
977 bhp_limit,
978 this->connections_.rho(),
979 this->getALQ(well_state),
980 thp_limit,
981 deferred_logger);
982 if ( (this->isProducer() && thp < thp_limit) || (this->isInjector() && thp > thp_limit) ) {
983 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
984 }
985 }
986 } else {
987 // defaulted BHP and there is a THP constraint
988 // default BHP limit is about 1 atm.
989 // when applied the hydrostatic pressure correction dp,
990 // most likely we get a negative value (bhp + dp)to search in the VFP table,
991 // which is not desirable.
992 // we assume we can operate under defaulted BHP limit and will violate the THP limit
993 // when operating under defaulted BHP limit.
994 this->operability_status_.operable_under_only_bhp_limit = true;
995 this->operability_status_.obey_thp_limit_under_bhp_limit = false;
996 }
997 }
998
999
1000
1001
1002
1003 template<typename TypeTag>
1004 void
1007 const WellState<Scalar>& well_state,
1008 DeferredLogger& deferred_logger)
1009 {
1010 const auto& summaryState = simulator.vanguard().summaryState();
1011 const auto obtain_bhp = this->isProducer() ? computeBhpAtThpLimitProd(well_state, simulator, summaryState, deferred_logger)
1012 : computeBhpAtThpLimitInj(simulator, summaryState, deferred_logger);
1013
1014 if (obtain_bhp) {
1015 this->operability_status_.can_obtain_bhp_with_thp_limit = true;
1016
1017 const Scalar bhp_limit = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1018 this->operability_status_.obey_bhp_limit_with_thp_limit = this->isProducer() ?
1019 *obtain_bhp >= bhp_limit : *obtain_bhp <= bhp_limit ;
1020
1021 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1022 if (this->isProducer() && *obtain_bhp < thp_limit) {
1023 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1024 + " bars is SMALLER than thp limit "
1025 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1026 + " bars as a producer for well " + name();
1027 deferred_logger.debug(msg);
1028 }
1029 else if (this->isInjector() && *obtain_bhp > thp_limit) {
1030 const std::string msg = " obtained bhp " + std::to_string(unit::convert::to(*obtain_bhp, unit::barsa))
1031 + " bars is LARGER than thp limit "
1032 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1033 + " bars as a injector for well " + name();
1034 deferred_logger.debug(msg);
1035 }
1036 } else {
1037 this->operability_status_.can_obtain_bhp_with_thp_limit = false;
1038 this->operability_status_.obey_bhp_limit_with_thp_limit = false;
1039 if (!this->wellIsStopped()) {
1040 const Scalar thp_limit = this->getTHPConstraint(summaryState);
1041 deferred_logger.debug(" could not find bhp value at thp limit "
1042 + std::to_string(unit::convert::to(thp_limit, unit::barsa))
1043 + " bar for well " + name() + ", the well might need to be closed ");
1044 }
1045 }
1046 }
1047
1048
1049
1050
1051
1052 template<typename TypeTag>
1053 bool
1055 allDrawDownWrongDirection(const Simulator& simulator) const
1056 {
1057 bool all_drawdown_wrong_direction = true;
1058
1059 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
1060 const int cell_idx = this->well_cells_[perf];
1061 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/0);
1062 const auto& fs = intQuants.fluidState();
1063
1064 const Scalar pressure = this->getPerfCellPressure(fs).value();
1065 const Scalar bhp = this->primary_variables_.eval(Bhp).value();
1066
1067 // Pressure drawdown (also used to determine direction of flow)
1068 const Scalar well_pressure = bhp + this->connections_.pressure_diff(perf);
1069 const Scalar drawdown = pressure - well_pressure;
1070
1071 // for now, if there is one perforation can produce/inject in the correct
1072 // direction, we consider this well can still produce/inject.
1073 // TODO: it can be more complicated than this to cause wrong-signed rates
1074 if ( (drawdown < 0. && this->isInjector()) ||
1075 (drawdown > 0. && this->isProducer()) ) {
1076 all_drawdown_wrong_direction = false;
1077 break;
1078 }
1079 }
1080
1081 const auto& comm = this->parallel_well_info_.communication();
1082 if (comm.size() > 1)
1083 {
1084 all_drawdown_wrong_direction =
1085 (comm.min(all_drawdown_wrong_direction ? 1 : 0) == 1);
1086 }
1087
1088 return all_drawdown_wrong_direction;
1089 }
1090
1091
1092
1093
1094 template<typename TypeTag>
1095 bool
1098 const WellState<Scalar>& well_state,
1099 DeferredLogger& deferred_logger)
1100 {
1101 const Scalar bhp = well_state.well(this->index_of_well_).bhp;
1102 std::vector<Scalar> well_rates;
1103 computeWellRatesWithBhp(simulator, bhp, well_rates, deferred_logger);
1104
1105 const Scalar sign = (this->isProducer()) ? -1. : 1.;
1106 const Scalar threshold = sign * std::numeric_limits<Scalar>::min();
1107
1108 bool can_produce_inject = false;
1109 for (const auto value : well_rates) {
1110 if (this->isProducer() && value < threshold) {
1111 can_produce_inject = true;
1112 break;
1113 } else if (this->isInjector() && value > threshold) {
1114 can_produce_inject = true;
1115 break;
1116 }
1117 }
1118
1119 if (!can_produce_inject) {
1120 deferred_logger.debug(" well " + name() + " CANNOT produce or inejct ");
1121 }
1122
1123 return can_produce_inject;
1124 }
1125
1126
1127
1128
1129
1130 template<typename TypeTag>
1131 bool
1133 openCrossFlowAvoidSingularity(const Simulator& simulator) const
1134 {
1135 return !this->getAllowCrossFlow() && allDrawDownWrongDirection(simulator);
1136 }
1137
1138
1139
1140
1141 template<typename TypeTag>
1145 const WellState<Scalar>& well_state) const
1146 {
1147 auto prop_func = typename StdWellEval::StdWellConnections::PressurePropertyFunctions {
1148 // getTemperature
1149 [&model = simulator.model()](int cell_idx, int phase_idx)
1150 {
1151 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1152 .fluidState().temperature(phase_idx).value();
1153 },
1154
1155 // getSaltConcentration
1156 [&model = simulator.model()](int cell_idx)
1157 {
1158 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1159 .fluidState().saltConcentration().value();
1160 },
1161
1162 // getPvtRegionIdx
1163 [&model = simulator.model()](int cell_idx)
1164 {
1165 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1166 .fluidState().pvtRegionIndex();
1167 }
1168 };
1169
1170 if constexpr (Indices::enableSolvent) {
1171 prop_func.solventInverseFormationVolumeFactor =
1172 [&model = simulator.model()](int cell_idx)
1173 {
1174 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1175 .solventInverseFormationVolumeFactor().value();
1176 };
1177
1178 prop_func.solventRefDensity = [&model = simulator.model()](int cell_idx)
1179 {
1180 return model.intensiveQuantities(cell_idx, /* time_idx = */ 0)
1181 .solventRefDensity();
1182 };
1183 }
1184
1185 return this->connections_.computePropertiesForPressures(well_state, prop_func);
1186 }
1187
1188
1189
1190
1191
1192 template<typename TypeTag>
1195 getWellConvergence(const Simulator& simulator,
1196 const WellState<Scalar>& well_state,
1197 const std::vector<Scalar>& B_avg,
1198 DeferredLogger& deferred_logger,
1199 const bool relax_tolerance) const
1200 {
1201 // the following implementation assume that the polymer is always after the w-o-g phases
1202 // For the polymer, energy and foam cases, there is one more mass balance equations of reservoir than wells
1203 assert((int(B_avg.size()) == this->num_components_) || has_polymer || has_energy || has_foam || has_brine || has_zFraction || has_micp);
1204
1205 Scalar tol_wells = this->param_.tolerance_wells_;
1206 // use stricter tolerance for stopped wells and wells under zero rate target control.
1207 constexpr Scalar stopped_factor = 1.e-4;
1208 // use stricter tolerance for dynamic thp to ameliorate network convergence
1209 constexpr Scalar dynamic_thp_factor = 1.e-1;
1210 if (this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger)) {
1211 tol_wells = tol_wells*stopped_factor;
1212 } else if (this->getDynamicThpLimit()) {
1213 tol_wells = tol_wells*dynamic_thp_factor;
1214 }
1215
1216 std::vector<Scalar> res;
1217 ConvergenceReport report = this->StdWellEval::getWellConvergence(well_state,
1218 B_avg,
1219 this->param_.max_residual_allowed_,
1220 tol_wells,
1221 this->param_.relaxed_tolerance_flow_well_,
1222 relax_tolerance,
1223 this->wellIsStopped(),
1224 res,
1225 deferred_logger);
1226
1227 checkConvergenceExtraEqs(res, report);
1228
1229 return report;
1230 }
1231
1232
1233
1234
1235
1236 template<typename TypeTag>
1237 void
1239 updateProductivityIndex(const Simulator& simulator,
1240 const WellProdIndexCalculator<Scalar>& wellPICalc,
1241 WellState<Scalar>& well_state,
1242 DeferredLogger& deferred_logger) const
1243 {
1244 auto fluidState = [&simulator, this](const int perf)
1245 {
1246 const auto cell_idx = this->well_cells_[perf];
1247 return simulator.model()
1248 .intensiveQuantities(cell_idx, /*timeIdx=*/ 0).fluidState();
1249 };
1250
1251 const int np = this->number_of_phases_;
1252 auto setToZero = [np](Scalar* x) -> void
1253 {
1254 std::fill_n(x, np, 0.0);
1255 };
1256
1257 auto addVector = [np](const Scalar* src, Scalar* dest) -> void
1258 {
1259 std::transform(src, src + np, dest, dest, std::plus<>{});
1260 };
1261
1262 auto& ws = well_state.well(this->index_of_well_);
1263 auto& perf_data = ws.perf_data;
1264 auto* wellPI = ws.productivity_index.data();
1265 auto* connPI = perf_data.prod_index.data();
1266
1267 setToZero(wellPI);
1268
1269 const auto preferred_phase = this->well_ecl_.getPreferredPhase();
1270 auto subsetPerfID = 0;
1271
1272 for (const auto& perf : *this->perf_data_) {
1273 auto allPerfID = perf.ecl_index;
1274
1275 auto connPICalc = [&wellPICalc, allPerfID](const Scalar mobility) -> Scalar
1276 {
1277 return wellPICalc.connectionProdIndStandard(allPerfID, mobility);
1278 };
1279
1280 std::vector<Scalar> mob(this->num_components_, 0.0);
1281 getMobility(simulator, static_cast<int>(subsetPerfID), mob, deferred_logger);
1282
1283 const auto& fs = fluidState(subsetPerfID);
1284 setToZero(connPI);
1285
1286 if (this->isInjector()) {
1287 this->computeConnLevelInjInd(fs, preferred_phase, connPICalc,
1288 mob, connPI, deferred_logger);
1289 }
1290 else { // Production or zero flow rate
1291 this->computeConnLevelProdInd(fs, connPICalc, mob, connPI);
1292 }
1293
1294 addVector(connPI, wellPI);
1295
1296 ++subsetPerfID;
1297 connPI += np;
1298 }
1299
1300 // Sum with communication in case of distributed well.
1301 const auto& comm = this->parallel_well_info_.communication();
1302 if (comm.size() > 1) {
1303 comm.sum(wellPI, np);
1304 }
1305
1306 assert ((static_cast<int>(subsetPerfID) == this->number_of_perforations_) &&
1307 "Internal logic error in processing connections for PI/II");
1308 }
1309
1310
1311
1312 template<typename TypeTag>
1315 const WellState<Scalar>& well_state,
1316 const WellConnectionProps& props,
1317 DeferredLogger& deferred_logger)
1318 {
1319 // Cell level dynamic property call-back functions as fall-back
1320 // option for calculating connection level mixture densities in
1321 // stopped or zero-rate producer wells.
1322 const auto prop_func = typename StdWellEval::StdWellConnections::DensityPropertyFunctions {
1323 // This becomes slightly more palatable with C++20's designated
1324 // initialisers.
1325
1326 // mobility: Phase mobilities in specified cell.
1327 [&model = simulator.model()](const int cell,
1328 const std::vector<int>& phases,
1329 std::vector<Scalar>& mob)
1330 {
1331 const auto& iq = model.intensiveQuantities(cell, /* time_idx = */ 0);
1332
1333 std::transform(phases.begin(), phases.end(), mob.begin(),
1334 [&iq](const int phase) { return iq.mobility(phase).value(); });
1335 },
1336
1337 // densityInCell: Reservoir condition phase densities in
1338 // specified cell.
1339 [&model = simulator.model()](const int cell,
1340 const std::vector<int>& phases,
1341 std::vector<Scalar>& rho)
1342 {
1343 const auto& fs = model.intensiveQuantities(cell, /* time_idx = */ 0).fluidState();
1344
1345 std::transform(phases.begin(), phases.end(), rho.begin(),
1346 [&fs](const int phase) { return fs.density(phase).value(); });
1347 }
1348 };
1349
1350 const auto stopped_or_zero_rate_target = this->
1351 stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1352
1353 this->connections_
1354 .computeProperties(stopped_or_zero_rate_target, well_state,
1355 prop_func, props, deferred_logger);
1356 }
1357
1358
1359
1360
1361
1362 template<typename TypeTag>
1363 void
1366 const WellState<Scalar>& well_state,
1367 DeferredLogger& deferred_logger)
1368 {
1369 const auto props = computePropertiesForWellConnectionPressures
1370 (simulator, well_state);
1371
1372 computeWellConnectionDensitesPressures(simulator, well_state,
1373 props, deferred_logger);
1374 }
1375
1376
1377
1378
1379
1380 template<typename TypeTag>
1381 void
1383 solveEqAndUpdateWellState(const Simulator& simulator,
1384 WellState<Scalar>& well_state,
1385 DeferredLogger& deferred_logger)
1386 {
1387 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1388
1389 // We assemble the well equations, then we check the convergence,
1390 // which is why we do not put the assembleWellEq here.
1391 BVectorWell dx_well(1);
1392 dx_well[0].resize(this->primary_variables_.numWellEq());
1393 this->linSys_.solve( dx_well);
1394
1395 updateWellState(simulator, dx_well, well_state, deferred_logger);
1396 }
1397
1398
1399
1400
1401
1402 template<typename TypeTag>
1403 void
1406 const WellState<Scalar>& well_state,
1407 DeferredLogger& deferred_logger)
1408 {
1409 updatePrimaryVariables(simulator, well_state, deferred_logger);
1410 initPrimaryVariablesEvaluation();
1411 computeWellConnectionPressures(simulator, well_state, deferred_logger);
1412 this->computeAccumWell();
1413 }
1414
1415
1416
1417 template<typename TypeTag>
1418 void
1420 apply(const BVector& x, BVector& Ax) const
1421 {
1422 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1423
1424 if (this->param_.matrix_add_well_contributions_)
1425 {
1426 // Contributions are already in the matrix itself
1427 return;
1428 }
1429
1430 this->linSys_.apply(x, Ax);
1431 }
1432
1433
1434
1435
1436 template<typename TypeTag>
1437 void
1439 apply(BVector& r) const
1440 {
1441 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1442
1443 this->linSys_.apply(r);
1444 }
1445
1446
1447
1448
1449 template<typename TypeTag>
1450 void
1453 const BVector& x,
1454 WellState<Scalar>& well_state,
1455 DeferredLogger& deferred_logger)
1456 {
1457 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1458
1459 BVectorWell xw(1);
1460 xw[0].resize(this->primary_variables_.numWellEq());
1461
1462 this->linSys_.recoverSolutionWell(x, xw);
1463 updateWellState(simulator, xw, well_state, deferred_logger);
1464 }
1465
1466
1467
1468
1469 template<typename TypeTag>
1470 void
1472 computeWellRatesWithBhp(const Simulator& simulator,
1473 const Scalar& bhp,
1474 std::vector<Scalar>& well_flux,
1475 DeferredLogger& deferred_logger) const
1476 {
1477
1478 const int np = this->number_of_phases_;
1479 well_flux.resize(np, 0.0);
1480
1481 const bool allow_cf = this->getAllowCrossFlow();
1482
1483 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
1484 const int cell_idx = this->well_cells_[perf];
1485 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1486 // flux for each perforation
1487 std::vector<Scalar> mob(this->num_components_, 0.);
1488 getMobility(simulator, perf, mob, deferred_logger);
1489 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
1490 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1491 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
1492
1493 std::vector<Scalar> cq_s(this->num_components_, 0.);
1494 PerforationRates<Scalar> perf_rates;
1495 computePerfRate(intQuants, mob, bhp, Tw, perf, allow_cf,
1496 cq_s, perf_rates, deferred_logger);
1497
1498 for(int p = 0; p < np; ++p) {
1499 well_flux[this->modelCompIdxToFlowCompIdx(p)] += cq_s[p];
1500 }
1501
1502 // the solvent contribution is added to the gas potentials
1503 if constexpr (has_solvent) {
1504 const auto& pu = this->phaseUsage();
1505 assert(pu.phase_used[Gas]);
1506 const int gas_pos = pu.phase_pos[Gas];
1507 well_flux[gas_pos] += cq_s[Indices::contiSolventEqIdx];
1508 }
1509 }
1510 this->parallel_well_info_.communication().sum(well_flux.data(), well_flux.size());
1511 }
1512
1513
1514
1515 template<typename TypeTag>
1516 void
1519 const Scalar& bhp,
1520 std::vector<Scalar>& well_flux,
1521 DeferredLogger& deferred_logger) const
1522 {
1523 // creating a copy of the well itself, to avoid messing up the explicit information
1524 // during this copy, the only information not copied properly is the well controls
1525 StandardWell<TypeTag> well_copy(*this);
1526
1527 // iterate to get a more accurate well density
1528 // create a copy of the well_state to use. If the operability checking is sucessful, we use this one
1529 // to replace the original one
1530 WellState<Scalar> well_state_copy = simulator.problem().wellModel().wellState();
1531 const auto& group_state = simulator.problem().wellModel().groupState();
1532
1533 // Get the current controls.
1534 const auto& summary_state = simulator.vanguard().summaryState();
1535 auto inj_controls = well_copy.well_ecl_.isInjector()
1536 ? well_copy.well_ecl_.injectionControls(summary_state)
1537 : Well::InjectionControls(0);
1538 auto prod_controls = well_copy.well_ecl_.isProducer()
1539 ? well_copy.well_ecl_.productionControls(summary_state) :
1540 Well::ProductionControls(0);
1541
1542 // Set current control to bhp, and bhp value in state, modify bhp limit in control object.
1543 auto& ws = well_state_copy.well(this->index_of_well_);
1544 if (well_copy.well_ecl_.isInjector()) {
1545 inj_controls.bhp_limit = bhp;
1546 ws.injection_cmode = Well::InjectorCMode::BHP;
1547 } else {
1548 prod_controls.bhp_limit = bhp;
1549 ws.production_cmode = Well::ProducerCMode::BHP;
1550 }
1551 ws.bhp = bhp;
1552
1553 // initialized the well rates with the potentials i.e. the well rates based on bhp
1554 const int np = this->number_of_phases_;
1555 const Scalar sign = this->well_ecl_.isInjector() ? 1.0 : -1.0;
1556 for (int phase = 0; phase < np; ++phase){
1557 well_state_copy.wellRates(this->index_of_well_)[phase]
1558 = sign * ws.well_potentials[phase];
1559 }
1560 well_copy.updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
1562 well_copy.computeAccumWell();
1563
1564 const double dt = simulator.timeStepSize();
1565 const bool converged = well_copy.iterateWellEqWithControl(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1566 if (!converged) {
1567 const std::string msg = " well " + name() + " did not get converged during well potential calculations "
1568 " potentials are computed based on unconverged solution";
1569 deferred_logger.debug(msg);
1570 }
1571 well_copy.updatePrimaryVariables(simulator, well_state_copy, deferred_logger);
1572 well_copy.computeWellConnectionPressures(simulator, well_state_copy, deferred_logger);
1574 well_copy.computeWellRatesWithBhp(simulator, bhp, well_flux, deferred_logger);
1575 }
1576
1577
1578
1579
1580 template<typename TypeTag>
1581 std::vector<typename StandardWell<TypeTag>::Scalar>
1584 DeferredLogger& deferred_logger,
1585 const WellState<Scalar>& well_state) const
1586 {
1587 std::vector<Scalar> potentials(this->number_of_phases_, 0.0);
1588 const auto& summary_state = simulator.vanguard().summaryState();
1589
1590 const auto& well = this->well_ecl_;
1591 if (well.isInjector()){
1592 const auto& controls = this->well_ecl_.injectionControls(summary_state);
1593 auto bhp_at_thp_limit = computeBhpAtThpLimitInj(simulator, summary_state, deferred_logger);
1594 if (bhp_at_thp_limit) {
1595 const Scalar bhp = std::min(*bhp_at_thp_limit,
1596 static_cast<Scalar>(controls.bhp_limit));
1597 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1598 } else {
1599 deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL",
1600 "Failed in getting converged thp based potential calculation for well "
1601 + name() + ". Instead the bhp based value is used");
1602 const Scalar bhp = controls.bhp_limit;
1603 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1604 }
1605 } else {
1606 computeWellRatesWithThpAlqProd(
1607 simulator, summary_state,
1608 deferred_logger, potentials, this->getALQ(well_state)
1609 );
1610 }
1611
1612 return potentials;
1613 }
1614
1615 template<typename TypeTag>
1616 bool
1619 std::vector<Scalar>& well_potentials,
1620 DeferredLogger& deferred_logger) const
1621 {
1622 // Create a copy of the well.
1623 // TODO: check if we can avoid taking multiple copies. Call from updateWellPotentials
1624 // is allready a copy, but not from other calls.
1625 StandardWell<TypeTag> well_copy(*this);
1626
1627 // store a copy of the well state, we don't want to update the real well state
1628 WellState<Scalar> well_state_copy = simulator.problem().wellModel().wellState();
1629 const auto& group_state = simulator.problem().wellModel().groupState();
1630 auto& ws = well_state_copy.well(this->index_of_well_);
1631
1632 // get current controls
1633 const auto& summary_state = simulator.vanguard().summaryState();
1634 auto inj_controls = well_copy.well_ecl_.isInjector()
1635 ? well_copy.well_ecl_.injectionControls(summary_state)
1636 : Well::InjectionControls(0);
1637 auto prod_controls = well_copy.well_ecl_.isProducer()
1638 ? well_copy.well_ecl_.productionControls(summary_state) :
1639 Well::ProductionControls(0);
1640
1641 // prepare/modify well state and control
1642 well_copy.prepareForPotentialCalculations(summary_state, well_state_copy, inj_controls, prod_controls);
1643
1644 // update connection pressures relative to updated bhp to get better estimate of connection dp
1645 const int num_perf = ws.perf_data.size();
1646 for (int perf = 0; perf < num_perf; ++perf) {
1647 ws.perf_data.pressure[perf] = ws.bhp + well_copy.connections_.pressure_diff(perf);
1648 }
1649 // initialize rates from previous potentials
1650 const int np = this->number_of_phases_;
1651 bool trivial = true;
1652 for (int phase = 0; phase < np; ++phase){
1653 trivial = trivial && (ws.well_potentials[phase] == 0.0) ;
1654 }
1655 if (!trivial) {
1656 const Scalar sign = well_copy.well_ecl_.isInjector() ? 1.0 : -1.0;
1657 for (int phase = 0; phase < np; ++phase) {
1658 ws.surface_rates[phase] = sign * ws.well_potentials[phase];
1659 }
1660 }
1661
1662 well_copy.calculateExplicitQuantities(simulator, well_state_copy, deferred_logger);
1663 const double dt = simulator.timeStepSize();
1664 // iterate to get a solution at the given bhp.
1665 bool converged = false;
1666 if (this->well_ecl_.isProducer() && this->wellHasTHPConstraints(summary_state)) {
1667 converged = well_copy.solveWellWithTHPConstraint(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1668 } else {
1669 converged = well_copy.iterateWellEqWithSwitching(simulator, dt, inj_controls, prod_controls, well_state_copy, group_state, deferred_logger);
1670 }
1671
1672 // fetch potentials (sign is updated on the outside).
1673 well_potentials.clear();
1674 well_potentials.resize(np, 0.0);
1675 for (int comp_idx = 0; comp_idx < this->num_components_; ++comp_idx) {
1676 if (has_solvent && comp_idx == Indices::contiSolventEqIdx) continue; // we do not store the solvent in the well_potentials
1677 const EvalWell rate = well_copy.primary_variables_.getQs(comp_idx);
1678 well_potentials[this->modelCompIdxToFlowCompIdx(comp_idx)] = rate.value();
1679 }
1680
1681 // the solvent contribution is added to the gas potentials
1682 if constexpr (has_solvent) {
1683 const auto& pu = this->phaseUsage();
1684 assert(pu.phase_used[Gas]);
1685 const int gas_pos = pu.phase_pos[Gas];
1686 const EvalWell rate = well_copy.primary_variables_.getQs(Indices::contiSolventEqIdx);
1687 well_potentials[gas_pos] += rate.value();
1688 }
1689 return converged;
1690 }
1691
1692
1693 template<typename TypeTag>
1697 const SummaryState &summary_state,
1698 DeferredLogger& deferred_logger,
1699 std::vector<Scalar>& potentials,
1700 Scalar alq) const
1701 {
1702 Scalar bhp;
1703 auto bhp_at_thp_limit = computeBhpAtThpLimitProdWithAlq(
1704 simulator, summary_state, alq, deferred_logger);
1705 if (bhp_at_thp_limit) {
1706 const auto& controls = this->well_ecl_.productionControls(summary_state);
1707 bhp = std::max(*bhp_at_thp_limit,
1708 static_cast<Scalar>(controls.bhp_limit));
1709 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1710 }
1711 else {
1712 deferred_logger.warning("FAILURE_GETTING_CONVERGED_POTENTIAL",
1713 "Failed in getting converged thp based potential calculation for well "
1714 + name() + ". Instead the bhp based value is used");
1715 const auto& controls = this->well_ecl_.productionControls(summary_state);
1716 bhp = controls.bhp_limit;
1717 computeWellRatesWithBhp(simulator, bhp, potentials, deferred_logger);
1718 }
1719 return bhp;
1720 }
1721
1722 template<typename TypeTag>
1723 void
1726 const SummaryState& summary_state,
1727 DeferredLogger& deferred_logger,
1728 std::vector<Scalar>& potentials,
1729 Scalar alq) const
1730 {
1731 /*double bhp =*/
1732 computeWellRatesAndBhpWithThpAlqProd(simulator,
1733 summary_state,
1734 deferred_logger,
1735 potentials,
1736 alq);
1737 }
1738
1739 template<typename TypeTag>
1740 void
1742 computeWellPotentials(const Simulator& simulator,
1743 const WellState<Scalar>& well_state,
1744 std::vector<Scalar>& well_potentials,
1745 DeferredLogger& deferred_logger) // const
1746 {
1747 const auto [compute_potential, bhp_controlled_well] =
1748 this->WellInterfaceGeneric<Scalar>::computeWellPotentials(well_potentials, well_state);
1749
1750 if (!compute_potential) {
1751 return;
1752 }
1753
1754 bool converged_implicit = false;
1755 if (this->param_.local_well_solver_control_switching_) {
1756 converged_implicit = computeWellPotentialsImplicit(simulator, well_potentials, deferred_logger);
1757 }
1758 if (!converged_implicit) {
1759 // does the well have a THP related constraint?
1760 const auto& summaryState = simulator.vanguard().summaryState();
1761 if (!Base::wellHasTHPConstraints(summaryState) || bhp_controlled_well) {
1762 // get the bhp value based on the bhp constraints
1763 Scalar bhp = WellBhpThpCalculator(*this).mostStrictBhpFromBhpLimits(summaryState);
1764
1765 // In some very special cases the bhp pressure target are
1766 // temporary violated. This may lead to too small or negative potentials
1767 // that could lead to premature shutting of wells.
1768 // As a remedy the bhp that gives the largest potential is used.
1769 // For converged cases, ws.bhp <=bhp for injectors and ws.bhp >= bhp,
1770 // and the potentials will be computed using the limit as expected.
1771 const auto& ws = well_state.well(this->index_of_well_);
1772 if (this->isInjector())
1773 bhp = std::max(ws.bhp, bhp);
1774 else
1775 bhp = std::min(ws.bhp, bhp);
1776
1777 assert(std::abs(bhp) != std::numeric_limits<Scalar>::max());
1778 computeWellRatesWithBhpIterations(simulator, bhp, well_potentials, deferred_logger);
1779 } else {
1780 // the well has a THP related constraint
1781 well_potentials = computeWellPotentialWithTHP(simulator, deferred_logger, well_state);
1782 }
1783 }
1784
1785 this->checkNegativeWellPotentials(well_potentials,
1786 this->param_.check_well_operability_,
1787 deferred_logger);
1788 }
1789
1790
1791
1792
1793
1794
1795
1796 template<typename TypeTag>
1799 connectionDensity([[maybe_unused]] const int globalConnIdx,
1800 const int openConnIdx) const
1801 {
1802 return (openConnIdx < 0)
1803 ? 0.0
1804 : this->connections_.rho(openConnIdx);
1805 }
1806
1807
1808
1809
1810
1811 template<typename TypeTag>
1812 void
1814 updatePrimaryVariables(const Simulator& simulator,
1815 const WellState<Scalar>& well_state,
1816 DeferredLogger& deferred_logger)
1817 {
1818 if (!this->isOperableAndSolvable() && !this->wellIsStopped()) return;
1819
1820 const bool stop_or_zero_rate_target = this->stoppedOrZeroRateTarget(simulator, well_state, deferred_logger);
1821 this->primary_variables_.update(well_state, stop_or_zero_rate_target, deferred_logger);
1822
1823 // other primary variables related to polymer injection
1824 if constexpr (Base::has_polymermw) {
1825 this->primary_variables_.updatePolyMW(well_state);
1826 }
1827
1828 this->primary_variables_.checkFinite(deferred_logger);
1829 }
1830
1831
1832
1833
1834 template<typename TypeTag>
1837 getRefDensity() const
1838 {
1839 return this->connections_.rho();
1840 }
1841
1842
1843
1844
1845 template<typename TypeTag>
1846 void
1849 const int perf,
1850 std::vector<EvalWell>& mob,
1851 DeferredLogger& deferred_logger) const
1852 {
1853 const int cell_idx = this->well_cells_[perf];
1854 const auto& int_quant = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
1855 const EvalWell polymer_concentration = this->extendEval(int_quant.polymerConcentration());
1856
1857 // TODO: not sure should based on the well type or injecting/producing peforations
1858 // it can be different for crossflow
1859 if (this->isInjector()) {
1860 // assume fully mixing within injecting wellbore
1861 const auto& visc_mult_table = PolymerModule::plyviscViscosityMultiplierTable(int_quant.pvtRegionIndex());
1862 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
1863 mob[waterCompIdx] /= (this->extendEval(int_quant.waterViscosityCorrection()) * visc_mult_table.eval(polymer_concentration, /*extrapolate=*/true) );
1864 }
1865
1866 if (PolymerModule::hasPlyshlog()) {
1867 // we do not calculate the shear effects for injection wells when they do not
1868 // inject polymer.
1869 if (this->isInjector() && this->wpolymer() == 0.) {
1870 return;
1871 }
1872 // compute the well water velocity with out shear effects.
1873 // TODO: do we need to turn on crossflow here?
1874 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
1875 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
1876
1877 std::vector<EvalWell> cq_s(this->num_components_, {this->primary_variables_.numWellEq() + Indices::numEq, 0.});
1878 PerforationRates<Scalar> perf_rates;
1879 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(int_quant, cell_idx);
1880 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
1881 const std::vector<Scalar> Tw = this->wellIndex(perf, int_quant, trans_mult, wellstate_nupcol);
1882 computePerfRate(int_quant, mob, bhp, Tw, perf, allow_cf, cq_s,
1883 perf_rates, deferred_logger);
1884 // TODO: make area a member
1885 const Scalar area = 2 * M_PI * this->perf_rep_radius_[perf] * this->perf_length_[perf];
1886 const auto& material_law_manager = simulator.problem().materialLawManager();
1887 const auto& scaled_drainage_info =
1888 material_law_manager->oilWaterScaledEpsInfoDrainage(cell_idx);
1889 const Scalar swcr = scaled_drainage_info.Swcr;
1890 const EvalWell poro = this->extendEval(int_quant.porosity());
1891 const EvalWell sw = this->extendEval(int_quant.fluidState().saturation(FluidSystem::waterPhaseIdx));
1892 // guard against zero porosity and no water
1893 const EvalWell denom = max( (area * poro * (sw - swcr)), 1e-12);
1894 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
1895 EvalWell water_velocity = cq_s[waterCompIdx] / denom * this->extendEval(int_quant.fluidState().invB(FluidSystem::waterPhaseIdx));
1896
1897 if (PolymerModule::hasShrate()) {
1898 // the equation for the water velocity conversion for the wells and reservoir are from different version
1899 // of implementation. It can be changed to be more consistent when possible.
1900 water_velocity *= PolymerModule::shrate( int_quant.pvtRegionIndex() ) / this->bore_diameters_[perf];
1901 }
1902 const EvalWell shear_factor = PolymerModule::computeShearFactor(polymer_concentration,
1903 int_quant.pvtRegionIndex(),
1904 water_velocity);
1905 // modify the mobility with the shear factor.
1906 mob[waterCompIdx] /= shear_factor;
1907 }
1908 }
1909
1910 template<typename TypeTag>
1911 void
1913 {
1914 this->linSys_.extract(jacobian);
1915 }
1916
1917
1918 template <typename TypeTag>
1919 void
1921 const BVector& weights,
1922 const int pressureVarIndex,
1923 const bool use_well_weights,
1924 const WellState<Scalar>& well_state) const
1925 {
1926 this->linSys_.extractCPRPressureMatrix(jacobian,
1927 weights,
1928 pressureVarIndex,
1929 use_well_weights,
1930 *this,
1931 Bhp,
1932 well_state);
1933 }
1934
1935
1936
1937 template<typename TypeTag>
1940 pskinwater(const Scalar throughput,
1941 const EvalWell& water_velocity,
1942 DeferredLogger& deferred_logger) const
1943 {
1944 if constexpr (Base::has_polymermw) {
1945 const int water_table_id = this->polymerWaterTable_();
1946 if (water_table_id <= 0) {
1947 OPM_DEFLOG_THROW(std::runtime_error,
1948 fmt::format("Unused SKPRWAT table id used for well {}", name()),
1949 deferred_logger);
1950 }
1951 const auto& water_table_func = PolymerModule::getSkprwatTable(water_table_id);
1952 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
1953 // the skin pressure when injecting water, which also means the polymer concentration is zero
1954 EvalWell pskin_water(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
1955 pskin_water = water_table_func.eval(throughput_eval, water_velocity);
1956 return pskin_water;
1957 } else {
1958 OPM_DEFLOG_THROW(std::runtime_error,
1959 fmt::format("Polymermw is not activated, while injecting "
1960 "skin pressure is requested for well {}", name()),
1961 deferred_logger);
1962 }
1963 }
1964
1965
1966
1967
1968
1969 template<typename TypeTag>
1972 pskin(const Scalar throughput,
1973 const EvalWell& water_velocity,
1974 const EvalWell& poly_inj_conc,
1975 DeferredLogger& deferred_logger) const
1976 {
1977 if constexpr (Base::has_polymermw) {
1978 const Scalar sign = water_velocity >= 0. ? 1.0 : -1.0;
1979 const EvalWell water_velocity_abs = abs(water_velocity);
1980 if (poly_inj_conc == 0.) {
1981 return sign * pskinwater(throughput, water_velocity_abs, deferred_logger);
1982 }
1983 const int polymer_table_id = this->polymerTable_();
1984 if (polymer_table_id <= 0) {
1985 OPM_DEFLOG_THROW(std::runtime_error,
1986 fmt::format("Unavailable SKPRPOLY table id used for well {}", name()),
1987 deferred_logger);
1988 }
1989 const auto& skprpolytable = PolymerModule::getSkprpolyTable(polymer_table_id);
1990 const Scalar reference_concentration = skprpolytable.refConcentration;
1991 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
1992 // the skin pressure when injecting water, which also means the polymer concentration is zero
1993 EvalWell pskin_poly(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
1994 pskin_poly = skprpolytable.table_func.eval(throughput_eval, water_velocity_abs);
1995 if (poly_inj_conc == reference_concentration) {
1996 return sign * pskin_poly;
1997 }
1998 // poly_inj_conc != reference concentration of the table, then some interpolation will be required
1999 const EvalWell pskin_water = pskinwater(throughput, water_velocity_abs, deferred_logger);
2000 const EvalWell pskin = pskin_water + (pskin_poly - pskin_water) / reference_concentration * poly_inj_conc;
2001 return sign * pskin;
2002 } else {
2003 OPM_DEFLOG_THROW(std::runtime_error,
2004 fmt::format("Polymermw is not activated, while injecting "
2005 "skin pressure is requested for well {}", name()),
2006 deferred_logger);
2007 }
2008 }
2009
2010
2011
2012
2013
2014 template<typename TypeTag>
2017 wpolymermw(const Scalar throughput,
2018 const EvalWell& water_velocity,
2019 DeferredLogger& deferred_logger) const
2020 {
2021 if constexpr (Base::has_polymermw) {
2022 const int table_id = this->polymerInjTable_();
2023 const auto& table_func = PolymerModule::getPlymwinjTable(table_id);
2024 const EvalWell throughput_eval(this->primary_variables_.numWellEq() + Indices::numEq, throughput);
2025 EvalWell molecular_weight(this->primary_variables_.numWellEq() + Indices::numEq, 0.);
2026 if (this->wpolymer() == 0.) { // not injecting polymer
2027 return molecular_weight;
2028 }
2029 molecular_weight = table_func.eval(throughput_eval, abs(water_velocity));
2030 return molecular_weight;
2031 } else {
2032 OPM_DEFLOG_THROW(std::runtime_error,
2033 fmt::format("Polymermw is not activated, while injecting "
2034 "polymer molecular weight is requested for well {}", name()),
2035 deferred_logger);
2036 }
2037 }
2038
2039
2040
2041
2042
2043 template<typename TypeTag>
2044 void
2046 updateWaterThroughput([[maybe_unused]] const double dt,
2047 WellState<Scalar>& well_state) const
2048 {
2049 if constexpr (Base::has_polymermw) {
2050 if (!this->isInjector()) {
2051 return;
2052 }
2053
2054 auto& perf_water_throughput = well_state.well(this->index_of_well_)
2055 .perf_data.water_throughput;
2056
2057 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
2058 const Scalar perf_water_vel =
2059 this->primary_variables_.value(Bhp + 1 + perf);
2060
2061 // we do not consider the formation damage due to water
2062 // flowing from reservoir into wellbore
2063 if (perf_water_vel > Scalar{0}) {
2064 perf_water_throughput[perf] += perf_water_vel * dt;
2065 }
2066 }
2067 }
2068 }
2069
2070
2071
2072
2073
2074 template<typename TypeTag>
2075 void
2077 handleInjectivityRate(const Simulator& simulator,
2078 const int perf,
2079 std::vector<EvalWell>& cq_s) const
2080 {
2081 const int cell_idx = this->well_cells_[perf];
2082 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2083 const auto& fs = int_quants.fluidState();
2084 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2085 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2086 const int wat_vel_index = Bhp + 1 + perf;
2087 const unsigned water_comp_idx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
2088
2089 // water rate is update to use the form from water velocity, since water velocity is
2090 // a primary variable now
2091 cq_s[water_comp_idx] = area * this->primary_variables_.eval(wat_vel_index) * b_w;
2092 }
2093
2094
2095
2096
2097 template<typename TypeTag>
2098 void
2100 handleInjectivityEquations(const Simulator& simulator,
2101 const WellState<Scalar>& well_state,
2102 const int perf,
2103 const EvalWell& water_flux_s,
2104 DeferredLogger& deferred_logger)
2105 {
2106 const int cell_idx = this->well_cells_[perf];
2107 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2108 const auto& fs = int_quants.fluidState();
2109 const EvalWell b_w = this->extendEval(fs.invB(FluidSystem::waterPhaseIdx));
2110 const EvalWell water_flux_r = water_flux_s / b_w;
2111 const Scalar area = M_PI * this->bore_diameters_[perf] * this->perf_length_[perf];
2112 const EvalWell water_velocity = water_flux_r / area;
2113 const int wat_vel_index = Bhp + 1 + perf;
2114
2115 // equation for the water velocity
2116 const EvalWell eq_wat_vel = this->primary_variables_.eval(wat_vel_index) - water_velocity;
2117
2118 const auto& ws = well_state.well(this->index_of_well_);
2119 const auto& perf_data = ws.perf_data;
2120 const auto& perf_water_throughput = perf_data.water_throughput;
2121 const Scalar throughput = perf_water_throughput[perf];
2122 const int pskin_index = Bhp + 1 + this->number_of_perforations_ + perf;
2123
2124 EvalWell poly_conc(this->primary_variables_.numWellEq() + Indices::numEq, 0.0);
2125 poly_conc.setValue(this->wpolymer());
2126
2127 // equation for the skin pressure
2128 const EvalWell eq_pskin = this->primary_variables_.eval(pskin_index)
2129 - pskin(throughput, this->primary_variables_.eval(wat_vel_index), poly_conc, deferred_logger);
2130
2132 assembleInjectivityEq(eq_pskin,
2133 eq_wat_vel,
2134 pskin_index,
2135 wat_vel_index,
2136 cell_idx,
2137 this->primary_variables_.numWellEq(),
2138 this->linSys_);
2139 }
2140
2141
2142
2143
2144
2145 template<typename TypeTag>
2146 void
2148 checkConvergenceExtraEqs(const std::vector<Scalar>& res,
2149 ConvergenceReport& report) const
2150 {
2151 // if different types of extra equations are involved, this function needs to be refactored further
2152
2153 // checking the convergence of the extra equations related to polymer injectivity
2154 if constexpr (Base::has_polymermw) {
2155 WellConvergence(*this).
2156 checkConvergencePolyMW(res, Bhp, this->param_.max_residual_allowed_, report);
2157 }
2158 }
2159
2160
2161
2162
2163
2164 template<typename TypeTag>
2165 void
2167 updateConnectionRatePolyMW(const EvalWell& cq_s_poly,
2168 const IntensiveQuantities& int_quants,
2169 const WellState<Scalar>& well_state,
2170 const int perf,
2171 std::vector<RateVector>& connectionRates,
2172 DeferredLogger& deferred_logger) const
2173 {
2174 // the source term related to transport of molecular weight
2175 EvalWell cq_s_polymw = cq_s_poly;
2176 if (this->isInjector()) {
2177 const int wat_vel_index = Bhp + 1 + perf;
2178 const EvalWell water_velocity = this->primary_variables_.eval(wat_vel_index);
2179 if (water_velocity > 0.) { // injecting
2180 const auto& ws = well_state.well(this->index_of_well_);
2181 const auto& perf_water_throughput = ws.perf_data.water_throughput;
2182 const Scalar throughput = perf_water_throughput[perf];
2183 const EvalWell molecular_weight = wpolymermw(throughput, water_velocity, deferred_logger);
2184 cq_s_polymw *= molecular_weight;
2185 } else {
2186 // we do not consider the molecular weight from the polymer
2187 // going-back to the wellbore through injector
2188 cq_s_polymw *= 0.;
2189 }
2190 } else if (this->isProducer()) {
2191 if (cq_s_polymw < 0.) {
2192 cq_s_polymw *= this->extendEval(int_quants.polymerMoleWeight() );
2193 } else {
2194 // we do not consider the molecular weight from the polymer
2195 // re-injecting back through producer
2196 cq_s_polymw *= 0.;
2197 }
2198 }
2199 connectionRates[perf][Indices::contiPolymerMWEqIdx] = Base::restrictEval(cq_s_polymw);
2200 }
2201
2202
2203
2204
2205
2206
2207 template<typename TypeTag>
2208 std::optional<typename StandardWell<TypeTag>::Scalar>
2211 const Simulator& simulator,
2212 const SummaryState& summary_state,
2213 DeferredLogger& deferred_logger) const
2214 {
2215 return computeBhpAtThpLimitProdWithAlq(simulator,
2216 summary_state,
2217 this->getALQ(well_state),
2218 deferred_logger);
2219 }
2220
2221 template<typename TypeTag>
2222 std::optional<typename StandardWell<TypeTag>::Scalar>
2225 const SummaryState& summary_state,
2226 const Scalar alq_value,
2227 DeferredLogger& deferred_logger) const
2228 {
2229 // Make the frates() function.
2230 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2231 // Not solving the well equations here, which means we are
2232 // calculating at the current Fg/Fw values of the
2233 // well. This does not matter unless the well is
2234 // crossflowing, and then it is likely still a good
2235 // approximation.
2236 std::vector<Scalar> rates(3);
2237 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2238 this->adaptRatesForVFP(rates);
2239 return rates;
2240 };
2241
2242 Scalar max_pressure = 0.0;
2243 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
2244 const int cell_idx = this->well_cells_[perf];
2245 const auto& int_quants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2246 const auto& fs = int_quants.fluidState();
2247 Scalar pressure_cell = this->getPerfCellPressure(fs).value();
2248 max_pressure = std::max(max_pressure, pressure_cell);
2249 }
2250 auto bhpAtLimit = WellBhpThpCalculator(*this).computeBhpAtThpLimitProd(frates,
2251 summary_state,
2252 max_pressure,
2253 this->connections_.rho(),
2254 alq_value,
2255 this->getTHPConstraint(summary_state),
2256 deferred_logger);
2257
2258 if (bhpAtLimit) {
2259 auto v = frates(*bhpAtLimit);
2260 if (std::all_of(v.cbegin(), v.cend(), [](Scalar i){ return i <= 0; }) ) {
2261 return bhpAtLimit;
2262 }
2263 }
2264
2265
2266 auto fratesIter = [this, &simulator, &deferred_logger](const Scalar bhp) {
2267 // Solver the well iterations to see if we are
2268 // able to get a solution with an update
2269 // solution
2270 std::vector<Scalar> rates(3);
2271 computeWellRatesWithBhpIterations(simulator, bhp, rates, deferred_logger);
2272 this->adaptRatesForVFP(rates);
2273 return rates;
2274 };
2275
2276 bhpAtLimit = WellBhpThpCalculator(*this).computeBhpAtThpLimitProd(fratesIter,
2277 summary_state,
2278 max_pressure,
2279 this->connections_.rho(),
2280 alq_value,
2281 this->getTHPConstraint(summary_state),
2282 deferred_logger);
2283
2284
2285 if (bhpAtLimit) {
2286 // should we use fratesIter here since fratesIter is used in computeBhpAtThpLimitProd above?
2287 auto v = frates(*bhpAtLimit);
2288 if (std::all_of(v.cbegin(), v.cend(), [](Scalar i){ return i <= 0; }) ) {
2289 return bhpAtLimit;
2290 }
2291 }
2292
2293 // we still don't get a valied solution.
2294 return std::nullopt;
2295 }
2296
2297
2298
2299 template<typename TypeTag>
2300 std::optional<typename StandardWell<TypeTag>::Scalar>
2302 computeBhpAtThpLimitInj(const Simulator& simulator,
2303 const SummaryState& summary_state,
2304 DeferredLogger& deferred_logger) const
2305 {
2306 // Make the frates() function.
2307 auto frates = [this, &simulator, &deferred_logger](const Scalar bhp) {
2308 // Not solving the well equations here, which means we are
2309 // calculating at the current Fg/Fw values of the
2310 // well. This does not matter unless the well is
2311 // crossflowing, and then it is likely still a good
2312 // approximation.
2313 std::vector<Scalar> rates(3);
2314 computeWellRatesWithBhp(simulator, bhp, rates, deferred_logger);
2315 return rates;
2316 };
2317
2318 return WellBhpThpCalculator(*this).computeBhpAtThpLimitInj(frates,
2319 summary_state,
2320 this->connections_.rho(),
2321 1e-6,
2322 50,
2323 true,
2324 deferred_logger);
2325 }
2326
2327
2328
2329
2330
2331 template<typename TypeTag>
2332 bool
2334 iterateWellEqWithControl(const Simulator& simulator,
2335 const double dt,
2336 const Well::InjectionControls& inj_controls,
2337 const Well::ProductionControls& prod_controls,
2338 WellState<Scalar>& well_state,
2339 const GroupState<Scalar>& group_state,
2340 DeferredLogger& deferred_logger)
2341 {
2342 const int max_iter = this->param_.max_inner_iter_wells_;
2343 int it = 0;
2344 bool converged;
2345 bool relax_convergence = false;
2346 this->regularize_ = false;
2347 do {
2348 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
2349
2350 if (it > this->param_.strict_inner_iter_wells_) {
2351 relax_convergence = true;
2352 this->regularize_ = true;
2353 }
2354
2355 auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
2356
2357 converged = report.converged();
2358 if (converged) {
2359 break;
2360 }
2361
2362 ++it;
2363 solveEqAndUpdateWellState(simulator, well_state, deferred_logger);
2364
2365 // TODO: when this function is used for well testing purposes, will need to check the controls, so that we will obtain convergence
2366 // under the most restrictive control. Based on this converged results, we can check whether to re-open the well. Either we refactor
2367 // this function or we use different functions for the well testing purposes.
2368 // We don't allow for switching well controls while computing well potentials and testing wells
2369 // updateWellControl(simulator, well_state, deferred_logger);
2370 initPrimaryVariablesEvaluation();
2371 } while (it < max_iter);
2372
2373 return converged;
2374 }
2375
2376
2377 template<typename TypeTag>
2378 bool
2380 iterateWellEqWithSwitching(const Simulator& simulator,
2381 const double dt,
2382 const Well::InjectionControls& inj_controls,
2383 const Well::ProductionControls& prod_controls,
2384 WellState<Scalar>& well_state,
2385 const GroupState<Scalar>& group_state,
2386 DeferredLogger& deferred_logger,
2387 const bool fixed_control /*false*/,
2388 const bool fixed_status /*false*/)
2389 {
2390 const int max_iter = this->param_.max_inner_iter_wells_;
2391 int it = 0;
2392 bool converged = false;
2393 bool relax_convergence = false;
2394 this->regularize_ = false;
2395 const auto& summary_state = simulator.vanguard().summaryState();
2396
2397 // Always take a few (more than one) iterations after a switch before allowing a new switch
2398 // The optimal number here is subject to further investigation, but it has been observerved
2399 // that unless this number is >1, we may get stuck in a cycle
2400 constexpr int min_its_after_switch = 4;
2401 int its_since_last_switch = min_its_after_switch;
2402 int switch_count= 0;
2403 // if we fail to solve eqs, we reset status/operability before leaving
2404 const auto well_status_orig = this->wellStatus_;
2405 const auto operability_orig = this->operability_status_;
2406 auto well_status_cur = well_status_orig;
2407 int status_switch_count = 0;
2408 // don't allow opening wells that are stopped from schedule or has a stopped well state
2409 const bool allow_open = this->well_ecl_.getStatus() == WellStatus::OPEN &&
2410 well_state.well(this->index_of_well_).status == WellStatus::OPEN;
2411 // don't allow switcing for wells under zero rate target or requested fixed status and control
2412 const bool allow_switching =
2413 !this->wellUnderZeroRateTarget(simulator, well_state, deferred_logger) &&
2414 (!fixed_control || !fixed_status) && allow_open;
2415
2416 bool changed = false;
2417 bool final_check = false;
2418 // well needs to be set operable or else solving/updating of re-opened wells is skipped
2419 this->operability_status_.resetOperability();
2420 this->operability_status_.solvable = true;
2421 do {
2422 its_since_last_switch++;
2423 if (allow_switching && its_since_last_switch >= min_its_after_switch){
2424 const Scalar wqTotal = this->primary_variables_.eval(WQTotal).value();
2425 changed = this->updateWellControlAndStatusLocalIteration(simulator, well_state, group_state,
2426 inj_controls, prod_controls, wqTotal,
2427 deferred_logger, fixed_control, fixed_status);
2428 if (changed){
2429 its_since_last_switch = 0;
2430 switch_count++;
2431 if (well_status_cur != this->wellStatus_) {
2432 well_status_cur = this->wellStatus_;
2433 status_switch_count++;
2434 }
2435 }
2436 if (!changed && final_check) {
2437 break;
2438 } else {
2439 final_check = false;
2440 }
2441 }
2442
2443 assembleWellEqWithoutIteration(simulator, dt, inj_controls, prod_controls, well_state, group_state, deferred_logger);
2444
2445 if (it > this->param_.strict_inner_iter_wells_) {
2446 relax_convergence = true;
2447 this->regularize_ = true;
2448 }
2449
2450 auto report = getWellConvergence(simulator, well_state, Base::B_avg_, deferred_logger, relax_convergence);
2451
2452 converged = report.converged();
2453 if (converged) {
2454 // if equations are sufficiently linear they might converge in less than min_its_after_switch
2455 // in this case, make sure all constraints are satisfied before returning
2456 if (switch_count > 0 && its_since_last_switch < min_its_after_switch) {
2457 final_check = true;
2458 its_since_last_switch = min_its_after_switch;
2459 } else {
2460 break;
2461 }
2462 }
2463
2464 ++it;
2465 solveEqAndUpdateWellState(simulator, well_state, deferred_logger);
2466 initPrimaryVariablesEvaluation();
2467
2468 } while (it < max_iter);
2469
2470 if (converged) {
2471 if (allow_switching){
2472 // update operability if status change
2473 const bool is_stopped = this->wellIsStopped();
2474 if (this->wellHasTHPConstraints(summary_state)){
2475 this->operability_status_.can_obtain_bhp_with_thp_limit = !is_stopped;
2476 this->operability_status_.obey_thp_limit_under_bhp_limit = !is_stopped;
2477 } else {
2478 this->operability_status_.operable_under_only_bhp_limit = !is_stopped;
2479 }
2480 }
2481 } else {
2482 this->wellStatus_ = well_status_orig;
2483 this->operability_status_ = operability_orig;
2484 const std::string message = fmt::format(" Well {} did not converge in {} inner iterations ("
2485 "{} switches, {} status changes).", this->name(), it, switch_count, status_switch_count);
2486 deferred_logger.debug(message);
2487 // add operability here as well ?
2488 }
2489 return converged;
2490 }
2491
2492 template<typename TypeTag>
2493 std::vector<typename StandardWell<TypeTag>::Scalar>
2495 computeCurrentWellRates(const Simulator& simulator,
2496 DeferredLogger& deferred_logger) const
2497 {
2498 // Calculate the rates that follow from the current primary variables.
2499 std::vector<Scalar> well_q_s(this->num_components_, 0.);
2500 const EvalWell& bhp = this->primary_variables_.eval(Bhp);
2501 const bool allow_cf = this->getAllowCrossFlow() || openCrossFlowAvoidSingularity(simulator);
2502 for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
2503 const int cell_idx = this->well_cells_[perf];
2504 const auto& intQuants = simulator.model().intensiveQuantities(cell_idx, /*timeIdx=*/ 0);
2505 std::vector<Scalar> mob(this->num_components_, 0.);
2506 getMobility(simulator, perf, mob, deferred_logger);
2507 std::vector<Scalar> cq_s(this->num_components_, 0.);
2508 Scalar trans_mult = simulator.problem().template wellTransMultiplier<Scalar>(intQuants, cell_idx);
2509 const auto& wellstate_nupcol = simulator.problem().wellModel().nupcolWellState().well(this->index_of_well_);
2510 const std::vector<Scalar> Tw = this->wellIndex(perf, intQuants, trans_mult, wellstate_nupcol);
2511 PerforationRates<Scalar> perf_rates;
2512 computePerfRate(intQuants, mob, bhp.value(), Tw, perf, allow_cf,
2513 cq_s, perf_rates, deferred_logger);
2514 for (int comp = 0; comp < this->num_components_; ++comp) {
2515 well_q_s[comp] += cq_s[comp];
2516 }
2517 }
2518 const auto& comm = this->parallel_well_info_.communication();
2519 if (comm.size() > 1)
2520 {
2521 comm.sum(well_q_s.data(), well_q_s.size());
2522 }
2523 return well_q_s;
2524 }
2525
2526
2527
2528 template <typename TypeTag>
2529 std::vector<typename StandardWell<TypeTag>::Scalar>
2531 getPrimaryVars() const
2532 {
2533 const int num_pri_vars = this->primary_variables_.numWellEq();
2534 std::vector<Scalar> retval(num_pri_vars);
2535 for (int ii = 0; ii < num_pri_vars; ++ii) {
2536 retval[ii] = this->primary_variables_.value(ii);
2537 }
2538 return retval;
2539 }
2540
2541
2542
2543
2544
2545 template <typename TypeTag>
2546 int
2548 setPrimaryVars(typename std::vector<Scalar>::const_iterator it)
2549 {
2550 const int num_pri_vars = this->primary_variables_.numWellEq();
2551 for (int ii = 0; ii < num_pri_vars; ++ii) {
2552 this->primary_variables_.setValue(ii, it[ii]);
2553 }
2554 return num_pri_vars;
2555 }
2556
2557
2558 template <typename TypeTag>
2561 connectionRateEnergy(const Scalar maxOilSaturation,
2562 const std::vector<EvalWell>& cq_s,
2563 const IntensiveQuantities& intQuants,
2564 DeferredLogger& deferred_logger) const
2565 {
2566 auto fs = intQuants.fluidState();
2567 Eval result = 0;
2568 for (unsigned phaseIdx = 0; phaseIdx < FluidSystem::numPhases; ++phaseIdx) {
2569 if (!FluidSystem::phaseIsActive(phaseIdx)) {
2570 continue;
2571 }
2572
2573 // convert to reservoir conditions
2574 EvalWell cq_r_thermal(this->primary_variables_.numWellEq() + Indices::numEq, 0.);
2575 const unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
2576 const bool both_oil_gas = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx) && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
2577 if (!both_oil_gas || FluidSystem::waterPhaseIdx == phaseIdx) {
2578 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2579 } else {
2580 // remove dissolved gas and vapporized oil
2581 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
2582 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
2583 // q_os = q_or * b_o + rv * q_gr * b_g
2584 // q_gs = q_gr * g_g + rs * q_or * b_o
2585 // q_gr = 1 / (b_g * d) * (q_gs - rs * q_os)
2586 // d = 1.0 - rs * rv
2587 const EvalWell d = this->extendEval(1.0 - fs.Rv() * fs.Rs());
2588 if (d <= 0.0) {
2589 deferred_logger.debug(
2590 fmt::format("Problematic d value {} obtained for well {}"
2591 " during calculateSinglePerf with rs {}"
2592 ", rv {}. Continue as if no dissolution (rs = 0) and"
2593 " vaporization (rv = 0) for this connection.",
2594 d, this->name(), fs.Rs(), fs.Rv()));
2595 cq_r_thermal = cq_s[activeCompIdx] / this->extendEval(fs.invB(phaseIdx));
2596 } else {
2597 if (FluidSystem::gasPhaseIdx == phaseIdx) {
2598 cq_r_thermal = (cq_s[gasCompIdx] -
2599 this->extendEval(fs.Rs()) * cq_s[oilCompIdx]) /
2600 (d * this->extendEval(fs.invB(phaseIdx)) );
2601 } else if (FluidSystem::oilPhaseIdx == phaseIdx) {
2602 // q_or = 1 / (b_o * d) * (q_os - rv * q_gs)
2603 cq_r_thermal = (cq_s[oilCompIdx] - this->extendEval(fs.Rv()) *
2604 cq_s[gasCompIdx]) /
2605 (d * this->extendEval(fs.invB(phaseIdx)) );
2606 }
2607 }
2608 }
2609
2610 // change temperature for injecting fluids
2611 if (this->isInjector() && cq_s[activeCompIdx] > 0.0){
2612 // only handles single phase injection now
2613 assert(this->well_ecl_.injectorType() != InjectorType::MULTI);
2614 fs.setTemperature(this->well_ecl_.inj_temperature());
2615 typedef typename std::decay<decltype(fs)>::type::Scalar FsScalar;
2616 typename FluidSystem::template ParameterCache<FsScalar> paramCache;
2617 const unsigned pvtRegionIdx = intQuants.pvtRegionIndex();
2618 paramCache.setRegionIndex(pvtRegionIdx);
2619 paramCache.setMaxOilSat(maxOilSaturation);
2620 paramCache.updatePhase(fs, phaseIdx);
2621
2622 const auto& rho = FluidSystem::density(fs, paramCache, phaseIdx);
2623 fs.setDensity(phaseIdx, rho);
2624 const auto& h = FluidSystem::enthalpy(fs, paramCache, phaseIdx);
2625 fs.setEnthalpy(phaseIdx, h);
2626 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2627 result += getValue(cq_r_thermal);
2628 } else {
2629 // compute the thermal flux
2630 cq_r_thermal *= this->extendEval(fs.enthalpy(phaseIdx)) * this->extendEval(fs.density(phaseIdx));
2631 result += Base::restrictEval(cq_r_thermal);
2632 }
2633 }
2634
2635 return result * this->well_efficiency_factor_;
2636 }
2637
2638
2639 template <typename TypeTag>
2640 template<class Value>
2641 void
2642 StandardWell<TypeTag>::
2643 gasOilPerfRateInj(const std::vector<Value>& cq_s,
2644 PerforationRates<Scalar>& perf_rates,
2645 const Value& rv,
2646 const Value& rs,
2647 const Value& pressure,
2648 const Value& rvw,
2649 DeferredLogger& deferred_logger) const
2650 {
2651 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
2652 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
2653 // TODO: the formulations here remain to be tested with cases with strong crossflow through production wells
2654 // s means standard condition, r means reservoir condition
2655 // q_os = q_or * b_o + rv * q_gr * b_g
2656 // q_gs = q_gr * b_g + rs * q_or * b_o
2657 // d = 1.0 - rs * rv
2658 // q_or = 1 / (b_o * d) * (q_os - rv * q_gs)
2659 // q_gr = 1 / (b_g * d) * (q_gs - rs * q_os)
2660
2661 const Scalar d = 1.0 - getValue(rv) * getValue(rs);
2662
2663 if (d <= 0.0) {
2664 deferred_logger.debug(dValueError(d, this->name(),
2665 "gasOilPerfRateInj",
2666 rs, rv, pressure));
2667 } else {
2668 // vaporized oil into gas
2669 // rv * q_gr * b_g = rv * (q_gs - rs * q_os) / d
2670 perf_rates.vap_oil = getValue(rv) * (getValue(cq_s[gasCompIdx]) - getValue(rs) * getValue(cq_s[oilCompIdx])) / d;
2671 // dissolved of gas in oil
2672 // rs * q_or * b_o = rs * (q_os - rv * q_gs) / d
2673 perf_rates.dis_gas = getValue(rs) * (getValue(cq_s[oilCompIdx]) - getValue(rv) * getValue(cq_s[gasCompIdx])) / d;
2674 }
2675
2676 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
2677 // q_ws = q_wr * b_w + rvw * q_gr * b_g
2678 // q_wr = 1 / b_w * (q_ws - rvw * q_gr * b_g) = 1 / b_w * (q_ws - rvw * 1 / d (q_gs - rs * q_os))
2679 // vaporized water in gas
2680 // rvw * q_gr * b_g = q_ws -q_wr *b_w = rvw * (q_gs -rs *q_os) / d
2681 perf_rates.vap_wat = getValue(rvw) * (getValue(cq_s[gasCompIdx]) - getValue(rs) * getValue(cq_s[oilCompIdx])) / d;
2682 }
2683 }
2684
2685
2686
2687 template <typename TypeTag>
2688 template<class Value>
2689 void
2690 StandardWell<TypeTag>::
2691 gasOilPerfRateProd(std::vector<Value>& cq_s,
2692 PerforationRates<Scalar>& perf_rates,
2693 const Value& rv,
2694 const Value& rs,
2695 const Value& rvw) const
2696 {
2697 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
2698 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
2699 const Value cq_sOil = cq_s[oilCompIdx];
2700 const Value cq_sGas = cq_s[gasCompIdx];
2701 const Value dis_gas = rs * cq_sOil;
2702 const Value vap_oil = rv * cq_sGas;
2703
2704 cq_s[gasCompIdx] += dis_gas;
2705 cq_s[oilCompIdx] += vap_oil;
2706
2707 // recording the perforation solution gas rate and solution oil rates
2708 if (this->isProducer()) {
2709 perf_rates.dis_gas = getValue(dis_gas);
2710 perf_rates.vap_oil = getValue(vap_oil);
2711 }
2712
2713 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
2714 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
2715 const Value vap_wat = rvw * cq_sGas;
2716 cq_s[waterCompIdx] += vap_wat;
2717 if (this->isProducer())
2718 perf_rates.vap_wat = getValue(vap_wat);
2719 }
2720 }
2721
2722
2723 template <typename TypeTag>
2724 template<class Value>
2725 void
2726 StandardWell<TypeTag>::
2727 gasWaterPerfRateProd(std::vector<Value>& cq_s,
2728 PerforationRates<Scalar>& perf_rates,
2729 const Value& rvw,
2730 const Value& rsw) const
2731 {
2732 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
2733 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
2734 const Value cq_sWat = cq_s[waterCompIdx];
2735 const Value cq_sGas = cq_s[gasCompIdx];
2736 const Value vap_wat = rvw * cq_sGas;
2737 const Value dis_gas_wat = rsw * cq_sWat;
2738 cq_s[waterCompIdx] += vap_wat;
2739 cq_s[gasCompIdx] += dis_gas_wat;
2740 if (this->isProducer()) {
2741 perf_rates.vap_wat = getValue(vap_wat);
2742 perf_rates.dis_gas_in_water = getValue(dis_gas_wat);
2743 }
2744 }
2745
2746
2747 template <typename TypeTag>
2748 template<class Value>
2749 void
2750 StandardWell<TypeTag>::
2751 gasWaterPerfRateInj(const std::vector<Value>& cq_s,
2752 PerforationRates<Scalar>& perf_rates,
2753 const Value& rvw,
2754 const Value& rsw,
2755 const Value& pressure,
2756 DeferredLogger& deferred_logger) const
2757
2758 {
2759 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
2760 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
2761
2762 const Scalar dw = 1.0 - getValue(rvw) * getValue(rsw);
2763
2764 if (dw <= 0.0) {
2765 deferred_logger.debug(dValueError(dw, this->name(),
2766 "gasWaterPerfRateInj",
2767 rsw, rvw, pressure));
2768 } else {
2769 // vaporized water into gas
2770 // rvw * q_gr * b_g = rvw * (q_gs - rsw * q_ws) / dw
2771 perf_rates.vap_wat = getValue(rvw) * (getValue(cq_s[gasCompIdx]) - getValue(rsw) * getValue(cq_s[waterCompIdx])) / dw;
2772 // dissolved gas in water
2773 // rsw * q_wr * b_w = rsw * (q_ws - rvw * q_gs) / dw
2774 perf_rates.dis_gas_in_water = getValue(rsw) * (getValue(cq_s[waterCompIdx]) - getValue(rvw) * getValue(cq_s[gasCompIdx])) / dw;
2775 }
2776 }
2777
2778
2779 template <typename TypeTag>
2780 template<class Value>
2781 void
2782 StandardWell<TypeTag>::
2783 disOilVapWatVolumeRatio(Value& volumeRatio,
2784 const Value& rvw,
2785 const Value& rsw,
2786 const Value& pressure,
2787 const std::vector<Value>& cmix_s,
2788 const std::vector<Value>& b_perfcells_dense,
2789 DeferredLogger& deferred_logger) const
2790 {
2791 const unsigned waterCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
2792 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
2793 // Incorporate RSW/RVW factors if both water and gas active
2794 const Value d = 1.0 - rvw * rsw;
2795
2796 if (d <= 0.0) {
2797 deferred_logger.debug(dValueError(d, this->name(),
2798 "disOilVapWatVolumeRatio",
2799 rsw, rvw, pressure));
2800 }
2801 const Value tmp_wat = d > 0.0 ? (cmix_s[waterCompIdx] - rvw * cmix_s[gasCompIdx]) / d
2802 : cmix_s[waterCompIdx];
2803 volumeRatio += tmp_wat / b_perfcells_dense[waterCompIdx];
2804
2805 const Value tmp_gas = d > 0.0 ? (cmix_s[gasCompIdx] - rsw * cmix_s[waterCompIdx]) / d
2806 : cmix_s[gasCompIdx];
2807 volumeRatio += tmp_gas / b_perfcells_dense[gasCompIdx];
2808 }
2809
2810
2811 template <typename TypeTag>
2812 template<class Value>
2813 void
2814 StandardWell<TypeTag>::
2815 gasOilVolumeRatio(Value& volumeRatio,
2816 const Value& rv,
2817 const Value& rs,
2818 const Value& pressure,
2819 const std::vector<Value>& cmix_s,
2820 const std::vector<Value>& b_perfcells_dense,
2821 DeferredLogger& deferred_logger) const
2822 {
2823 const unsigned oilCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::oilCompIdx);
2824 const unsigned gasCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::gasCompIdx);
2825 // Incorporate RS/RV factors if both oil and gas active
2826 const Value d = 1.0 - rv * rs;
2827
2828 if (d <= 0.0) {
2829 deferred_logger.debug(dValueError(d, this->name(),
2830 "gasOilVolumeRatio",
2831 rs, rv, pressure));
2832 }
2833 const Value tmp_oil = d > 0.0? (cmix_s[oilCompIdx] - rv * cmix_s[gasCompIdx]) / d : cmix_s[oilCompIdx];
2834 volumeRatio += tmp_oil / b_perfcells_dense[oilCompIdx];
2835
2836 const Value tmp_gas = d > 0.0? (cmix_s[gasCompIdx] - rs * cmix_s[oilCompIdx]) / d : cmix_s[gasCompIdx];
2837 volumeRatio += tmp_gas / b_perfcells_dense[gasCompIdx];
2838 }
2839
2840} // namespace Opm
#define OPM_DEFLOG_THROW(Exception, message, deferred_logger)
Definition: DeferredLoggingErrorHelpers.hpp:45
#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 debug(const std::string &tag, const std::string &message)
Definition: GroupState.hpp:38
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:186
Class handling assemble of the equation system for StandardWell.
Definition: StandardWellAssemble.hpp:43
Scalar pressure_diff(const unsigned perf) const
Returns pressure drop for a given perforation.
Definition: StandardWellConnections.hpp:105
StdWellConnections connections_
Connection level values.
Definition: StandardWellEval.hpp:107
PrimaryVariables primary_variables_
Primary variables for well.
Definition: StandardWellEval.hpp:101
Definition: StandardWell.hpp:59
void computeWellConnectionPressures(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1365
EvalWell wpolymermw(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2017
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: StandardWell_impl.hpp:1742
typename StdWellEval::EvalWell EvalWell
Definition: StandardWell.hpp:119
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &ebos_simulator, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2224
typename StdWellEval::BVectorWell BVectorWell
Definition: StandardWell.hpp:120
void checkOperabilityUnderBHPLimit(const WellState< Scalar > &well_state, const Simulator &simulator, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:937
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: StandardWell_impl.hpp:2334
void computeWellRatesWithThpAlqProd(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1725
void addWellContributions(SparseMatrixAdapter &mat) const override
Definition: StandardWell_impl.hpp:1912
std::vector< Scalar > getPrimaryVars() const override
Definition: StandardWell_impl.hpp:2531
void initPrimaryVariablesEvaluation() override
Definition: StandardWell_impl.hpp:113
void updateIPRImplicit(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:861
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1239
void updateWaterMobilityWithPolymer(const Simulator &simulator, const int perf, std::vector< EvalWell > &mob_water, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1848
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: StandardWell_impl.hpp:2380
StandardWell(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: StandardWell_impl.hpp:72
void updateWaterThroughput(const double dt, WellState< Scalar > &well_state) const override
Definition: StandardWell_impl.hpp:2046
std::vector< Scalar > computeCurrentWellRates(const Simulator &ebosSimulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:2495
virtual void init(const PhaseUsage *phase_usage_arg, const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const int num_cells, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step) override
Definition: StandardWell_impl.hpp:96
void updatePrimaryVariablesNewton(const BVectorWell &dwells, const bool stop_or_zero_rate_target, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:725
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
void calculateExplicitQuantities(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1405
void solveEqAndUpdateWellState(const Simulator &simulator, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1383
bool computeWellPotentialsImplicit(const Simulator &ebos_simulator, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1618
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1452
typename StdWellEval::StdWellConnections::Properties WellConnectionProps
Definition: StandardWell.hpp:270
void computeWellRatesWithBhp(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1472
WellConnectionProps computePropertiesForWellConnectionPressures(const Simulator &simulator, const WellState< Scalar > &well_state) const
Definition: StandardWell_impl.hpp:1144
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:647
void computeWellRatesWithBhpIterations(const Simulator &ebosSimulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:1518
void updateIPR(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: StandardWell_impl.hpp:768
std::vector< Scalar > computeWellPotentialWithTHP(const Simulator &ebosSimulator, DeferredLogger &deferred_logger, const WellState< Scalar > &well_state) const
Definition: StandardWell_impl.hpp:1583
void computePerfRate(const IntensiveQuantities &intQuants, const std::vector< Value > &mob, const Value &bhp, const std::vector< Scalar > &Tw, const int perf, const bool allow_cf, std::vector< Value > &cq_s, PerforationRates< Scalar > &perf_rates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:126
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:703
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: StandardWell_impl.hpp:344
void computeWellConnectionDensitesPressures(const Simulator &simulator, const WellState< Scalar > &well_state, const WellConnectionProps &props, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1314
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: StandardWell_impl.hpp:1420
std::optional< Scalar > computeBhpAtThpLimitProd(const WellState< Scalar > &well_state, const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2210
void checkConvergenceExtraEqs(const std::vector< Scalar > &res, ConvergenceReport &report) const
Definition: StandardWell_impl.hpp:2148
typename StdWellEval::Eval Eval
Definition: StandardWell.hpp:118
Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator &ebos_simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger, std::vector< Scalar > &potentials, Scalar alq) const
Definition: StandardWell_impl.hpp:1696
void updateConnectionRatePolyMW(const EvalWell &cq_s_poly, const IntensiveQuantities &int_quants, const WellState< Scalar > &well_state, const int perf, std::vector< RateVector > &connectionRates, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2167
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1133
void updatePrimaryVariables(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1814
virtual 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: StandardWell_impl.hpp:1195
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &simulator, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:2302
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: StandardWell_impl.hpp:1055
EvalWell pskin(const Scalar throughput, const EvalWell &water_velocity, const EvalWell &poly_inj_conc, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1972
void updateWellStateFromPrimaryVariables(const bool stop_or_zero_rate_target, WellState< Scalar > &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:748
static constexpr int numWellConservationEq
Definition: StandardWell.hpp:95
void handleInjectivityEquations(const Simulator &simulator, const WellState< Scalar > &well_state, const int perf, const EvalWell &water_flux_s, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:2100
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: StandardWell_impl.hpp:2548
void assembleWellEqWithoutIterationImpl(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: StandardWell_impl.hpp:370
void checkOperabilityUnderTHPLimit(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger) override
Definition: StandardWell_impl.hpp:1006
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellState< Scalar > &well_state) const override
Definition: StandardWell_impl.hpp:1920
EvalWell pskinwater(const Scalar throughput, const EvalWell &water_velocity, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:1940
bool canProduceInjectWithCurrentBhp(const Simulator &simulator, const WellState< Scalar > &well_state, DeferredLogger &deferred_logger)
Definition: StandardWell_impl.hpp:1097
void handleInjectivityRate(const Simulator &simulator, const int perf, std::vector< EvalWell > &cq_s) const
Definition: StandardWell_impl.hpp:2077
Scalar getRefDensity() const override
Definition: StandardWell_impl.hpp:1837
void calculateSinglePerf(const Simulator &simulator, const int perf, WellState< Scalar > &well_state, std::vector< RateVector > &connectionRates, std::vector< EvalWell > &cq_s, EvalWell &water_flux_s, EvalWell &cq_s_zfrac_effective, DeferredLogger &deferred_logger) const
Definition: StandardWell_impl.hpp:497
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: StandardWell_impl.hpp:1799
EvalWell getQs(const int compIdx) const
Returns scaled rate for a component.
Class for computing BHP limits.
Definition: WellBhpThpCalculator.hpp:41
std::optional< Scalar > computeBhpAtThpLimitProd(const std::function< std::vector< Scalar >(const Scalar)> &frates, const SummaryState &summary_state, const Scalar maxPerfPress, const Scalar rho, const Scalar alq_value, const Scalar thp_limit, DeferredLogger &deferred_logger) const
Compute BHP from THP limit for a producer.
std::optional< Scalar > computeBhpAtThpLimitInj(const std::function< std::vector< Scalar >(const Scalar)> &frates, const SummaryState &summary_state, const Scalar rho, const Scalar flo_rel_tol, const int max_iteration, const bool throwOnError, DeferredLogger &deferred_logger) const
Compute BHP from THP limit for an injector.
Scalar calculateThpFromBhp(const std::vector< Scalar > &rates, const Scalar bhp, const Scalar rho, const std::optional< Scalar > &alq, const Scalar thp_limit, DeferredLogger &deferred_logger) const
Calculates THP from BHP.
Scalar mostStrictBhpFromBhpLimits(const SummaryState &summaryState) const
Obtain the most strict BHP from BHP limits.
Definition: WellConvergence.hpp:38
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)
const int num_components_
Definition: WellInterfaceGeneric.hpp:282
Well well_ecl_
Definition: WellInterfaceGeneric.hpp:273
Definition: WellInterfaceIndices.hpp:34
Definition: WellInterface.hpp:73
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:78
void getMobility(const Simulator &simulator, const int perf, std::vector< Value > &mob, Callback &extendEval, DeferredLogger &deferred_logger) const
Definition: WellInterface_impl.hpp:1773
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:100
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:97
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:82
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:79
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:96
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:80
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:81
bool solveWellWithTHPConstraint(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:514
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:84
Definition: WellProdIndexCalculator.hpp:37
Scalar connectionProdIndStandard(const std::size_t connIdx, const Scalar connMobility) const
Definition: WellState.hpp:62
const SingleWellState< Scalar > & well(std::size_t well_index) const
Definition: WellState.hpp:307
std::vector< Scalar > & wellRates(std::size_t well_index)
One rate per well and phase.
Definition: WellState.hpp:272
@ NONE
Definition: DeferredLogger.hpp:46
Definition: blackoilboundaryratevector.hh:37
PhaseUsage phaseUsage(const Phases &phases)
Determine the active phases.
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:147
Static data associated with a well perforation.
Definition: PerforationData.hpp:30
Definition: PerforationData.hpp:41
Scalar dis_gas
Definition: PerforationData.hpp:42
Scalar vap_wat
Definition: PerforationData.hpp:45
Scalar vap_oil
Definition: PerforationData.hpp:44
Scalar dis_gas_in_water
Definition: PerforationData.hpp:43
Definition: BlackoilPhases.hpp:46