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
81 using Base::has_solvent;
82 using Base::has_zFraction;
83 using Base::has_polymer;
84 using Base::has_polymermw;
85 using Base::has_foam;
86 using Base::has_brine;
87 using Base::has_energy;
88 using Base::has_bioeffects;
89 using Base::has_micp;
90
91 using PolymerModule = BlackOilPolymerModule<TypeTag>;
92 using FoamModule = BlackOilFoamModule<TypeTag>;
93 using typename Base::PressureMatrix;
94
95 // number of the conservation equations
96 static constexpr int numWellConservationEq = Indices::numPhases + Indices::numSolvents;
97 // number of the well control equations
98 static constexpr int numWellControlEq = 1;
99 // number of the well equations that will always be used
100 // based on the solution strategy, there might be other well equations be introduced
101 static constexpr int numStaticWellEq = numWellConservationEq + numWellControlEq;
102
103 // the index for Bhp in primary variables and also the index of well control equation
104 // they both will be the last one in their respective system.
105 // TODO: we should have indices for the well equations and well primary variables separately
106 static constexpr int Bhp = numStaticWellEq - numWellControlEq;
107
108 using StdWellEval::WQTotal;
109
110 using typename Base::Scalar;
111
112 using Base::name;
113 using Base::Water;
114 using Base::Oil;
115 using Base::Gas;
116
117 using typename Base::BVector;
118
119 using Eval = typename StdWellEval::Eval;
120 using EvalWell = typename StdWellEval::EvalWell;
121 using BVectorWell = typename StdWellEval::BVectorWell;
122
123 using IndexTraits = typename FluidSystem::IndexTraitsType;
124 using WellStateType = WellState<Scalar, IndexTraits>;
125
126 StandardWell(const Well& well,
127 const ParallelWellInfo<Scalar>& pw_info,
128 const int time_step,
129 const ModelParameters& param,
130 const RateConverterType& rate_converter,
131 const int pvtRegionIdx,
132 const int num_conservation_quantities,
133 const int num_phases,
134 const int index_of_well,
135 const std::vector<PerforationData<Scalar>>& perf_data);
136
137 virtual void init(const std::vector<Scalar>& depth_arg,
138 const Scalar gravity_arg,
139 const std::vector<Scalar>& B_avg,
140 const bool changed_to_open_this_step) override;
141
143 virtual ConvergenceReport getWellConvergence(const Simulator& simulator,
144 const WellStateType& well_state,
145 const std::vector<Scalar>& B_avg,
146 DeferredLogger& deferred_logger,
147 const bool relax_tolerance) const override;
148
150 virtual void apply(const BVector& x, BVector& Ax) const override;
152 virtual void apply(BVector& r) const override;
153
156 void recoverWellSolutionAndUpdateWellState(const Simulator& simulator,
157 const BVector& x,
158 WellStateType& well_state,
159 DeferredLogger& deferred_logger) override;
160
162 void computeWellPotentials(const Simulator& simulator,
163 const WellStateType& well_state,
164 std::vector<Scalar>& well_potentials,
165 DeferredLogger& deferred_logger) /* const */ override;
166
167 void updatePrimaryVariables(const Simulator& simulator,
168 const WellStateType& well_state,
169 DeferredLogger& deferred_logger) override;
170
171 void solveEqAndUpdateWellState(const Simulator& simulator,
172 WellStateType& well_state,
173 DeferredLogger& deferred_logger) override;
174
175 void calculateExplicitQuantities(const Simulator& simulator,
176 const WellStateType& well_state,
177 DeferredLogger& deferred_logger) override; // should be const?
178
179 void updateProductivityIndex(const Simulator& simulator,
180 const WellProdIndexCalculator<Scalar>& wellPICalc,
181 WellStateType& well_state,
182 DeferredLogger& deferred_logger) const override;
183
184 Scalar connectionDensity(const int globalConnIdx,
185 const int openConnIdx) const override;
186
187 void addWellContributions(SparseMatrixAdapter& mat) const override;
188
189 void addWellPressureEquations(PressureMatrix& mat,
190 const BVector& x,
191 const int pressureVarIndex,
192 const bool use_well_weights,
193 const WellStateType& well_state) const override;
194
195 // iterate well equations with the specified control until converged
196 bool iterateWellEqWithControl(const Simulator& simulator,
197 const double dt,
198 const Well::InjectionControls& inj_controls,
199 const Well::ProductionControls& prod_controls,
200 WellStateType& well_state,
201 const GroupState<Scalar>& group_state,
202 DeferredLogger& deferred_logger) override;
203
204 // iterate well equations including control switching
205 bool iterateWellEqWithSwitching(const Simulator& simulator,
206 const double dt,
207 const Well::InjectionControls& inj_controls,
208 const Well::ProductionControls& prod_controls,
209 WellStateType& well_state,
210 const GroupState<Scalar>& group_state,
211 DeferredLogger& deferred_logger,
212 const bool fixed_control = false,
213 const bool fixed_status = false) override;
214
215 /* returns BHP */
216 Scalar computeWellRatesAndBhpWithThpAlqProd(const Simulator& ebos_simulator,
217 const SummaryState &summary_state,
218 DeferredLogger& deferred_logger,
219 std::vector<Scalar>& potentials,
220 Scalar alq) const;
221
222 void computeWellRatesWithThpAlqProd(const Simulator& ebos_simulator,
223 const SummaryState& summary_state,
224 DeferredLogger& deferred_logger,
225 std::vector<Scalar>& potentials,
226 Scalar alq) const;
227
228 std::optional<Scalar>
229 computeBhpAtThpLimitProdWithAlq(const Simulator& ebos_simulator,
230 const SummaryState& summary_state,
231 const Scalar alq_value,
232 DeferredLogger& deferred_logger,
233 bool iterate_if_no_solution) const override;
234
235 void updateIPRImplicit(const Simulator& simulator,
236 WellStateType& well_state,
237 DeferredLogger& deferred_logger) override;
238
239 void computeWellRatesWithBhp(const Simulator& ebosSimulator,
240 const Scalar& bhp,
241 std::vector<Scalar>& well_flux,
242 DeferredLogger& deferred_logger) const override;
243
244 // NOTE: These cannot be protected since they are used by GasLiftRuntime
245 using Base::vfp_properties_;
246
247 std::vector<Scalar>
248 computeCurrentWellRates(const Simulator& ebosSimulator,
249 DeferredLogger& deferred_logger) const override;
250
251 std::vector<Scalar> getPrimaryVars() const override;
252
253 int setPrimaryVars(typename std::vector<Scalar>::const_iterator it) override;
254
255 protected:
256 bool regularize_;
257
258 // updating the well_state based on well solution dwells
259 void updateWellState(const Simulator& simulator,
260 const BVectorWell& dwells,
261 WellStateType& well_state,
262 DeferredLogger& deferred_logger);
263
264 using WellConnectionProps = typename StdWellEval::StdWellConnections::Properties;
265
266 // Compute connection level PVT properties needed to calulate the
267 // pressure difference between well connections.
268 WellConnectionProps
269 computePropertiesForWellConnectionPressures(const Simulator& simulator,
270 const WellStateType& well_state) const;
271
272 void computeWellConnectionDensitesPressures(const Simulator& simulator,
273 const WellStateType& well_state,
274 const WellConnectionProps& props,
275 DeferredLogger& deferred_logger);
276
277 void computeWellConnectionPressures(const Simulator& simulator,
278 const WellStateType& well_state,
279 DeferredLogger& deferred_logger);
280
281 template<class Value>
282 void computePerfRate(const IntensiveQuantities& intQuants,
283 const std::vector<Value>& mob,
284 const Value& bhp,
285 const std::vector<Scalar>& Tw,
286 const int perf,
287 const bool allow_cf,
288 std::vector<Value>& cq_s,
289 PerforationRates<Scalar>& perf_rates,
290 DeferredLogger& deferred_logger) const;
291
292 template<class Value>
293 void computePerfRate(const std::vector<Value>& mob,
294 const Value& pressure,
295 const Value& bhp,
296 const Value& rs,
297 const Value& rv,
298 const Value& rvw,
299 const Value& rsw,
300 std::vector<Value>& b_perfcells_dense,
301 const std::vector<Scalar>& Tw,
302 const int perf,
303 const bool allow_cf,
304 const Value& skin_pressure,
305 const std::vector<Value>& cmix_s,
306 std::vector<Value>& cq_s,
307 PerforationRates<Scalar>& perf_rates,
308 DeferredLogger& deferred_logger) const;
309
310 void computeWellRatesWithBhpIterations(const Simulator& ebosSimulator,
311 const Scalar& bhp,
312 std::vector<Scalar>& well_flux,
313 DeferredLogger& deferred_logger) const override;
314
315 std::vector<Scalar>
316 computeWellPotentialWithTHP(const Simulator& ebosSimulator,
317 DeferredLogger& deferred_logger,
318 const WellStateType& well_state) const;
319
320 bool computeWellPotentialsImplicit(const Simulator& ebos_simulator,
321 const WellStateType& well_state,
322 std::vector<Scalar>& well_potentials,
323 DeferredLogger& deferred_logger) const;
324
325 // return the density at the perforation[0] of the rank owning this well,
326 // value is cached to minimize the number of broadcasts
327 Scalar getRefDensity() const override;
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 double dt,
351 const Well::InjectionControls& inj_controls,
352 const Well::ProductionControls& prod_controls,
353 WellStateType& well_state,
354 const GroupState<Scalar>& group_state,
355 DeferredLogger& deferred_logger) override;
356
357 void assembleWellEqWithoutIterationImpl(const Simulator& simulator,
358 const double dt,
359 const Well::InjectionControls& inj_controls,
360 const Well::ProductionControls& prod_controls,
361 WellStateType& well_state,
362 const GroupState<Scalar>& group_state,
363 DeferredLogger& deferred_logger);
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 DeferredLogger& deferred_logger) 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 SummaryState& summary_state,
449 DeferredLogger& deferred_logger) const;
450
451 std::optional<Scalar>
452 computeBhpAtThpLimitInj(const Simulator& simulator,
453 const SummaryState& summary_state,
454 DeferredLogger& deferred_logger) 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:58
Contains the high level supplements required to extend the black oil model by polymer.
Definition blackoilpolymermodules.hh:64
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
Definition GroupState.hpp:41
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:198
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:2507
virtual void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition StandardWell_impl.hpp:1415
virtual ConvergenceReport getWellConvergence(const Simulator &simulator, const WellStateType &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const bool relax_tolerance) const override
check whether the well equations get converged for this well
Definition StandardWell_impl.hpp:1186
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition StandardWell_impl.hpp:1737
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger) override
using the solution x to recover the solution xw for wells and applying xw to update Well State
Definition StandardWell_impl.hpp:1447
const std::string & name() const
Definition WellInterfaceGeneric.cpp:168
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:58
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:43
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