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