opm-simulators
Loading...
Searching...
No Matches
StandardWell.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
23#ifndef OPM_STANDARDWELL_HEADER_INCLUDED
24#define OPM_STANDARDWELL_HEADER_INCLUDED
25
26#include <opm/simulators/timestepping/ConvergenceReport.hpp>
28#include <opm/simulators/wells/RatioCalculator.hpp>
29#include <opm/simulators/wells/VFPInjProperties.hpp>
30#include <opm/simulators/wells/VFPProdProperties.hpp>
31#include <opm/simulators/wells/WellInterface.hpp>
32#include <opm/simulators/wells/WellProdIndexCalculator.hpp>
33#include <opm/simulators/wells/ParallelWellInfo.hpp>
34
41
42#include <opm/material/densead/Evaluation.hpp>
43#include <opm/input/eclipse/Schedule/ScheduleTypes.hpp>
44
45#include <opm/simulators/wells/StandardWellEval.hpp>
46
47#include <dune/common/dynvector.hh>
48#include <dune/common/dynmatrix.hh>
49
50#include <memory>
51#include <optional>
52
53namespace Opm
54{
55
56 template<typename TypeTag>
57 class StandardWell : public WellInterface<TypeTag>
58 , public StandardWellEval<GetPropType<TypeTag, Properties::FluidSystem>,
59 GetPropType<TypeTag, Properties::Indices>>
60 {
61
62 public:
63 using Base = WellInterface<TypeTag>;
64 using StdWellEval = StandardWellEval<GetPropType<TypeTag, Properties::FluidSystem>,
66
67 // TODO: some functions working with AD variables handles only with values (double) without
68 // dealing with derivatives. It can be beneficial to make functions can work with either AD or scalar value.
69 // And also, it can also be beneficial to make these functions hanle different types of AD variables.
70 using typename Base::Simulator;
71 using typename Base::IntensiveQuantities;
72 using typename Base::FluidSystem;
73 using typename Base::MaterialLaw;
74 using typename Base::ModelParameters;
75 using typename Base::Indices;
76 using typename Base::RateConverterType;
77 using typename Base::SparseMatrixAdapter;
78 using typename Base::FluidState;
79 using typename Base::RateVector;
80 using typename Base::GroupStateHelperType;
81
82 using Base::has_solvent;
83 using Base::has_zFraction;
84 using Base::has_polymer;
85 using Base::has_polymermw;
86 using Base::has_foam;
87 using Base::has_brine;
88 using Base::has_energy;
89 using Base::has_bioeffects;
90 using Base::has_micp;
91
92 using PolymerModule = BlackOilPolymerModule<TypeTag>;
93 using FoamModule = BlackOilFoamModule<TypeTag>;
94 using typename Base::PressureMatrix;
95
96 // number of the conservation equations
97 static constexpr int numWellConservationEq = Indices::numPhases + Indices::numSolvents;
98 // number of the well control equations
99 static constexpr int numWellControlEq = 1;
100 // number of the well equations that will always be used
101 // based on the solution strategy, there might be other well equations be introduced
102 static constexpr int numStaticWellEq = numWellConservationEq + numWellControlEq;
103
104 // the index for Bhp in primary variables and also the index of well control equation
105 // they both will be the last one in their respective system.
106 // TODO: we should have indices for the well equations and well primary variables separately
107 static constexpr int Bhp = numStaticWellEq - numWellControlEq;
108
109 using StdWellEval::WQTotal;
110
111 using typename Base::Scalar;
112
113 using Base::name;
114 using Base::Water;
115 using Base::Oil;
116 using Base::Gas;
117
118 using typename Base::BVector;
119
120 using Eval = typename StdWellEval::Eval;
121 using EvalWell = typename StdWellEval::EvalWell;
122 using BVectorWell = typename StdWellEval::BVectorWell;
123
124 using IndexTraits = typename FluidSystem::IndexTraitsType;
125 using WellStateType = WellState<Scalar, IndexTraits>;
126
127 StandardWell(const Well& well,
128 const ParallelWellInfo<Scalar>& pw_info,
129 const int time_step,
130 const ModelParameters& param,
131 const RateConverterType& rate_converter,
132 const int pvtRegionIdx,
133 const int num_conservation_quantities,
134 const int num_phases,
135 const int index_of_well,
136 const std::vector<PerforationData<Scalar>>& perf_data);
137
138 virtual void init(const std::vector<Scalar>& depth_arg,
139 const Scalar gravity_arg,
140 const std::vector<Scalar>& B_avg,
141 const bool changed_to_open_this_step) override;
142
144 virtual ConvergenceReport getWellConvergence(const GroupStateHelperType& groupStateHelper,
145 const std::vector<Scalar>& B_avg,
146 const bool relax_tolerance) const override;
147
149 virtual void apply(const BVector& x, BVector& Ax) const override;
151 virtual void apply(BVector& r) const override;
152
155 void recoverWellSolutionAndUpdateWellState(const Simulator& simulator,
156 const BVector& x,
157 const GroupStateHelperType& groupStateHelper,
158 WellStateType& well_state) override;
159
161 void computeWellPotentials(const Simulator& simulator,
162 const WellStateType& well_state,
163 const GroupStateHelperType& groupStateHelper,
164 std::vector<Scalar>& well_potentials) /* const */ override;
165
166 void updatePrimaryVariables(const GroupStateHelperType& groupStateHelper) override;
167
168 void solveEqAndUpdateWellState(const Simulator& simulator,
169 const GroupStateHelperType& groupStateHelper,
170 WellStateType& well_state) override;
171
172 void calculateExplicitQuantities(const Simulator& simulator,
173 const GroupStateHelperType& groupStateHelper) override; // should be const?
174
175 void updateProductivityIndex(const Simulator& simulator,
176 const WellProdIndexCalculator<Scalar>& wellPICalc,
177 WellStateType& well_state,
178 DeferredLogger& deferred_logger) const override;
179
180 Scalar connectionDensity(const int globalConnIdx,
181 const int openConnIdx) const override;
182
183 void addWellContributions(SparseMatrixAdapter& mat) const override;
184
185 void addWellPressureEquations(PressureMatrix& mat,
186 const BVector& x,
187 const int pressureVarIndex,
188 const bool use_well_weights,
189 const WellStateType& well_state) const override;
190
191 // iterate well equations with the specified control until converged
192 bool iterateWellEqWithControl(const Simulator& simulator,
193 const double dt,
194 const Well::InjectionControls& inj_controls,
195 const Well::ProductionControls& prod_controls,
196 const GroupStateHelperType& groupStateHelper,
197 WellStateType& well_state) override;
198
199 // iterate well equations including control switching
200 bool iterateWellEqWithSwitching(const Simulator& simulator,
201 const double dt,
202 const Well::InjectionControls& inj_controls,
203 const Well::ProductionControls& prod_controls,
204 const GroupStateHelperType& groupStateHelper,
205 WellStateType& well_state,
206 const bool fixed_control,
207 const bool fixed_status,
208 const bool solving_with_zero_rate) override;
209
210 /* returns BHP */
211 Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator& ebos_simulator,
212 const GroupStateHelperType& groupStateHelper,
213 const SummaryState &summary_state,
214 std::vector<Scalar>& potentials,
215 Scalar alq) const;
216
217 void computeWellRatesWithThpAlqProd(const Simulator& ebos_simulator,
218 const GroupStateHelperType& groupStateHelper,
219 const SummaryState& summary_state,
220 std::vector<Scalar>& potentials,
221 Scalar alq) const;
222
223 std::optional<Scalar>
224 computeBhpAtThpLimitProdWithAlq(const Simulator& ebos_simulator,
225 const GroupStateHelperType& groupStateHelper,
226 const SummaryState& summary_state,
227 const Scalar alq_value,
228 bool iterate_if_no_solution) const override;
229
230 void updateIPRImplicit(const Simulator& simulator,
231 const GroupStateHelperType& groupStateHelper,
232 WellStateType& well_state) override;
233
234 void computeWellRatesWithBhp(const Simulator& ebosSimulator,
235 const Scalar& bhp,
236 std::vector<Scalar>& well_flux,
237 DeferredLogger& deferred_logger) const override;
238
239 Scalar maxPerfPress(const Simulator& simulator) const override;
240
241 // NOTE: These cannot be protected since they are used by GasLiftRuntime
242 using Base::vfp_properties_;
243
244 std::vector<Scalar>
245 computeCurrentWellRates(const Simulator& ebosSimulator,
246 DeferredLogger& deferred_logger) const override;
247
248 std::vector<Scalar> getPrimaryVars() const override;
249
250 int setPrimaryVars(typename std::vector<Scalar>::const_iterator it) override;
251
252 protected:
253 bool regularize_;
254
255 // updating the well_state based on well solution dwells
256 void updateWellState(const Simulator& simulator,
257 const BVectorWell& dwells,
258 const GroupStateHelperType& groupStateHelper,
259 WellStateType& well_state);
260
261 using WellConnectionProps = typename StdWellEval::StdWellConnections::Properties;
262
263 // Compute connection level PVT properties needed to calulate the
264 // pressure difference between well connections.
265 WellConnectionProps
266 computePropertiesForWellConnectionPressures(const Simulator& simulator,
267 const WellStateType& well_state) const;
268
269 void computeWellConnectionDensitesPressures(const Simulator& simulator,
270 const GroupStateHelperType& groupStateHelper,
271 const WellConnectionProps& props);
272
273 void computeWellConnectionPressures(const Simulator& simulator,
274 const GroupStateHelperType& groupStateHelper);
275
276 template<class Value>
277 void computePerfRate(const IntensiveQuantities& intQuants,
278 const std::vector<Value>& mob,
279 const Value& bhp,
280 const std::vector<Value>& Tw,
281 const int perf,
282 const bool allow_cf,
283 std::vector<Value>& cq_s,
284 PerforationRates<Scalar>& perf_rates,
285 DeferredLogger& deferred_logger) const;
286
287 template<class Value>
288 void computePerfRate(const std::vector<Value>& mob,
289 const Value& pressure,
290 const Value& bhp,
291 const Value& rs,
292 const Value& rv,
293 const Value& rvw,
294 const Value& rsw,
295 std::vector<Value>& b_perfcells_dense,
296 const std::vector<Value>& Tw,
297 const int perf,
298 const bool allow_cf,
299 const Value& skin_pressure,
300 const std::vector<Value>& cmix_s,
301 std::vector<Value>& cq_s,
302 PerforationRates<Scalar>& perf_rates,
303 DeferredLogger& deferred_logger) const;
304
305 void computeWellRatesWithBhpIterations(const Simulator& ebosSimulator,
306 const Scalar& bhp,
307 const GroupStateHelperType& groupStateHelper,
308 std::vector<Scalar>& well_flux) const override;
309
310 std::vector<Scalar>
311 computeWellPotentialWithTHP(const Simulator& ebosSimulator,
312 const GroupStateHelperType& groupStateHelper,
313 const WellStateType& well_state) const;
314
315 bool computeWellPotentialsImplicit(const Simulator& ebos_simulator,
316 const GroupStateHelperType& groupStateHelper,
317 std::vector<Scalar>& well_potentials) const;
318
319 // return the density at the perforation[0] of the rank owning this well,
320 // value is cached to minimize the number of broadcasts
321 Scalar getRefDensity() const override;
322
323 // get the transmissibility multiplier for specific perforation
324 template<class Value>
325 void getTransMult(Value& trans_mult,
326 const Simulator& simulator,
327 const int cell_indx) const;
328
329 // get the mobility for specific perforation
330 template<class Value>
331 void getMobility(const Simulator& simulator,
332 const int perf,
333 std::vector<Value>& mob,
334 DeferredLogger& deferred_logger) const;
335
336 void updateWaterMobilityWithPolymer(const Simulator& simulator,
337 const int perf,
338 std::vector<EvalWell>& mob_water,
339 DeferredLogger& deferred_logger) const;
340
341 void updatePrimaryVariablesNewton(const BVectorWell& dwells,
342 const bool stop_or_zero_rate_target,
343 DeferredLogger& deferred_logger);
344
345 void updateWellStateFromPrimaryVariables(WellStateType& well_state,
346 const SummaryState& summary_state,
347 DeferredLogger& deferred_logger) const;
348
349 void assembleWellEqWithoutIteration(const Simulator& simulator,
350 const GroupStateHelperType& groupStateHelper,
351 const double dt,
352 const Well::InjectionControls& inj_controls,
353 const Well::ProductionControls& prod_controls,
354 WellStateType& well_state,
355 const bool solving_with_zero_rate) override;
356
357 void assembleWellEqWithoutIterationImpl(const Simulator& simulator,
358 const GroupStateHelperType& groupStateHelper,
359 const double dt,
360 const Well::InjectionControls& inj_controls,
361 const Well::ProductionControls& prod_controls,
362 WellStateType& well_state,
363 const bool solving_with_zero_rate);
364
365 void calculateSinglePerf(const Simulator& simulator,
366 const int perf,
367 WellStateType& well_state,
368 std::vector<RateVector>& connectionRates,
369 std::vector<EvalWell>& cq_s,
370 EvalWell& water_flux_s,
371 EvalWell& cq_s_zfrac_effective,
372 DeferredLogger& deferred_logger) const;
373
374 // check whether the well is operable under BHP limit with current reservoir condition
375 void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
376 const Simulator& simulator,
377 DeferredLogger& deferred_logger) override;
378
379 // check whether the well is operable under THP limit with current reservoir condition
380 void checkOperabilityUnderTHPLimit(const Simulator& simulator,
381 const WellStateType& well_state,
382 const GroupStateHelperType& groupStateHelper) override;
383
384 // updating the inflow based on the current reservoir condition
385 void updateIPR(const Simulator& simulator,
386 DeferredLogger& deferred_logger) const override;
387
388 // for a well, when all drawdown are in the wrong direction, then this well will not
389 // be able to produce/inject .
390 bool allDrawDownWrongDirection(const Simulator& simulator) const;
391
392 // turn on crossflow to avoid singular well equations
393 // when the well is banned from cross-flow and the BHP is not properly initialized,
394 // we turn on crossflow to avoid singular well equations. It can result in wrong-signed
395 // well rates, it can cause problem for THP calculation
396 // TODO: looking for better alternative to avoid wrong-signed well rates
397 bool openCrossFlowAvoidSingularity(const Simulator& simulator) const;
398
399 // calculate the skin pressure based on water velocity, throughput and polymer concentration.
400 // throughput is used to describe the formation damage during water/polymer injection.
401 // calculated skin pressure will be applied to the drawdown during perforation rate calculation
402 // to handle the effect from formation damage.
403 EvalWell pskin(const Scalar throughput,
404 const EvalWell& water_velocity,
405 const EvalWell& poly_inj_conc,
406 DeferredLogger& deferred_logger) const;
407
408 // calculate the skin pressure based on water velocity, throughput during water injection.
409 EvalWell pskinwater(const Scalar throughput,
410 const EvalWell& water_velocity,
411 DeferredLogger& deferred_logger) const;
412
413 // calculate the injecting polymer molecular weight based on the througput and water velocity
414 EvalWell wpolymermw(const Scalar throughput,
415 const EvalWell& water_velocity,
416 DeferredLogger& deferred_logger) const;
417
418 // modify the water rate for polymer injectivity study
419 void handleInjectivityRate(const Simulator& simulator,
420 const int perf,
421 std::vector<EvalWell>& cq_s) const;
422
423 // handle the extra equations for polymer injectivity study
424 void handleInjectivityEquations(const Simulator& simulator,
425 const WellStateType& well_state,
426 const int perf,
427 const EvalWell& water_flux_s,
428 DeferredLogger& deferred_logger);
429
430 void updateWaterThroughput(const double dt,
431 WellStateType& well_state) const override;
432
433 // checking convergence of extra equations, if there are any
434 void checkConvergenceExtraEqs(const std::vector<Scalar>& res,
435 ConvergenceReport& report) const;
436
437 // updating the connectionRates_ related polymer molecular weight
438 void updateConnectionRatePolyMW(const EvalWell& cq_s_poly,
439 const IntensiveQuantities& int_quants,
440 const WellStateType& well_state,
441 const int perf,
442 std::vector<RateVector>& connectionRates,
443 DeferredLogger& deferred_logger) const;
444
445 std::optional<Scalar>
446 computeBhpAtThpLimitProd(const WellStateType& well_state,
447 const Simulator& simulator,
448 const GroupStateHelperType& groupStateHelper,
449 const SummaryState& summary_state) const;
450
451 std::optional<Scalar>
452 computeBhpAtThpLimitInj(const Simulator& simulator,
453 const GroupStateHelperType& groupStateHelper,
454 const SummaryState& summary_state) const;
455
456 private:
457 Eval connectionRateEnergy(const std::vector<EvalWell>& cq_s,
458 const IntensiveQuantities& intQuants,
459 DeferredLogger& deferred_logger) const;
460
461 // density of the first perforation, might not be from this rank
462 Scalar cachedRefDensity{0};
463 };
464
465}
466
467#include "StandardWell_impl.hpp"
468
469#endif // OPM_STANDARDWELL_HEADER_INCLUDED
Facility for converting component rates at surface conditions to phase (voidage) rates at reservoir c...
Contains the classes required to extend the black-oil model by bioeffects.
Contains the classes required to extend the black-oil model by brine.
Contains the classes required to extend the black-oil model by solvent component.
Contains the classes required to extend the black-oil model to include the effects of foam.
Contains the classes required to extend the black-oil model by polymer.
Contains the classes required to extend the black-oil model by solvents.
Contains the high level supplements required to extend the black oil model to include the effects of ...
Definition blackoilfoammodules.hh:60
Contains the high level supplements required to extend the black oil model by polymer.
Definition blackoilpolymermodules.hh:65
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
virtual ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, const std::vector< Scalar > &B_avg, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition StandardWell_impl.hpp:1226
std::vector< Scalar > computeCurrentWellRates(const Simulator &ebosSimulator, DeferredLogger &deferred_logger) const override
Compute well rates based on current reservoir conditions and well variables.
Definition StandardWell_impl.hpp:2586
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials) override
computing the well potentials for group control
Definition StandardWell_impl.hpp:1796
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition StandardWell_impl.hpp:1453
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state) override
using the solution x to recover the solution xw for wells and applying xw to update Well State
Definition StandardWell_impl.hpp:1485
const std::string & name() const
Definition WellInterfaceGeneric.cpp:167
WellInterface(const Well &well, const ParallelWellInfo< Scalar > &pw_info, const int time_step, const ModelParameters &param, const RateConverterType &rate_converter, const int pvtRegionIdx, const int num_conservation_quantities, const int num_phases, const int index_of_well, const std::vector< PerforationData< Scalar > > &perf_data)
Constructor.
Definition WellInterface_impl.hpp:59
Collect per-connection static information to enable calculating connection-level or well-level produc...
Definition WellProdIndexCalculator.hpp:37
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:66
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:233
Static data associated with a well perforation.
Definition PerforationData.hpp:30
Definition PerforationData.hpp:41