BlackoilModelBase.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
3  Copyright 2014, 2015 Statoil ASA.
4  Copyright 2014, 2015 Dr. Markus Blatt - HPC-Simulation-Software & Services
5  Copyright 2015 NTNU
6 
7  This file is part of the Open Porous Media project (OPM).
8 
9  OPM is free software: you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation, either version 3 of the License, or
12  (at your option) any later version.
13 
14  OPM is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with OPM. If not, see <http://www.gnu.org/licenses/>.
21 */
22 
23 #ifndef OPM_BLACKOILMODELBASE_HEADER_INCLUDED
24 #define OPM_BLACKOILMODELBASE_HEADER_INCLUDED
25 
26 #include <cassert>
27 
35 #include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
36 
37 #include <array>
38 
39 struct Wells;
40 
41 namespace Opm {
42 
43  namespace parameter { class ParameterGroup; }
44  class DerivedGeology;
45  class RockCompressibility;
47  class VFPProperties;
48 
49 
52  {
54  explicit DefaultBlackoilSolutionState(const int np)
55  : pressure ( ADB::null())
56  , temperature( ADB::null())
57  , saturation(np, ADB::null())
58  , rs ( ADB::null())
59  , rv ( ADB::null())
60  , qs ( ADB::null())
61  , bhp ( ADB::null())
62  , canonical_phase_pressures(3, ADB::null())
63  {
64  }
65  ADB pressure;
67  std::vector<ADB> saturation;
68  ADB rs;
69  ADB rv;
70  ADB qs;
71  ADB bhp;
72  // Below are quantities stored in the state for optimization purposes.
73  std::vector<ADB> canonical_phase_pressures; // Always has 3 elements, even if only 2 phases active.
74  };
75 
76 
77 
78 
82  template <class ConcreteModel>
83  struct ModelTraits;
84 
85 
97  template<class Grid, class Implementation>
99  {
100  public:
101  // --------- Types and enums ---------
103  typedef ADB::V V;
104  typedef ADB::M M;
105 
110 
111  // --------- Public methods ---------
112 
128  BlackoilModelBase(const ModelParameters& param,
129  const Grid& grid ,
130  const BlackoilPropsAdInterface& fluid,
131  const DerivedGeology& geo ,
132  const RockCompressibility* rock_comp_props,
133  const Wells* wells,
134  const NewtonIterationBlackoilInterface& linsolver,
135  Opm::EclipseStateConstPtr eclState,
136  const bool has_disgas,
137  const bool has_vapoil,
138  const bool terminal_output);
139 
148  void setThresholdPressures(const std::vector<double>& threshold_pressures_by_face);
149 
154  void prepareStep(const double dt,
155  ReservoirState& reservoir_state,
156  WellState& well_state);
157 
163  void afterStep(const double dt,
164  ReservoirState& reservoir_state,
165  WellState& well_state);
166 
171  void assemble(const ReservoirState& reservoir_state,
172  WellState& well_state,
173  const bool initial_assembly);
174 
179  std::vector<double> computeResidualNorms() const;
180 
182  int sizeNonLinear() const;
183 
185  int linearIterationsLastSolve() const;
186 
189  V solveJacobianSystem() const;
190 
195  void updateState(const V& dx,
196  ReservoirState& reservoir_state,
197  WellState& well_state);
198 
200  bool terminalOutputEnabled() const;
201 
206  bool getConvergence(const double dt, const int iteration);
207 
209  int numPhases() const;
210 
214  int numMaterials() const;
215 
218  const std::string& materialName(int material_index) const;
219 
221  void updateEquationsScaling();
222 
223  protected:
224 
225  // --------- Types and enums ---------
226 
227  typedef Eigen::Array<double,
228  Eigen::Dynamic,
229  Eigen::Dynamic,
230  Eigen::RowMajor> DataBlock;
231 
234  std::vector<ADB> accum; // Accumulations
235  ADB mflux; // Mass flux (surface conditions)
236  ADB b; // Reciprocal FVF
237  ADB dh; // Pressure drop across int. interfaces
238  ADB mob; // Phase mobility (per cell)
239  };
240 
241  struct WellOps {
242  WellOps(const Wells* wells);
243  Eigen::SparseMatrix<double> w2p; // well -> perf (scatter)
244  Eigen::SparseMatrix<double> p2w; // perf -> well (gather)
245  };
246 
247  // --------- Data members ---------
248 
249  const Grid& grid_;
252  const RockCompressibility* rock_comp_props_;
253  const Wells* wells_;
256  // For each canonical phase -> true if active
257  const std::vector<bool> active_;
258  // Size = # active phases. Maps active -> canonical phase indices.
259  const std::vector<int> canph_;
260  const std::vector<int> cells_; // All grid cells
262  const WellOps wops_;
263  const bool has_disgas_;
264  const bool has_vapoil_;
265 
266  ModelParameters param_;
270 
271  std::vector<ReservoirResidualQuant> rq_;
272  std::vector<PhasePresence> phaseCondition_;
273  V isRs_;
274  V isRv_;
275  V isSg_;
276  V well_perforation_densities_; //Density of each well perforation
277  V well_perforation_pressure_diffs_; // Diff to bhp for each well perforation.
278 
280 
283 
284  std::vector<int> primalVariable_;
285  V pvdt_;
286  std::vector<std::string> material_name_;
287 
288  // --------- Protected methods ---------
289 
292  Implementation& asImpl()
293  {
294  return static_cast<Implementation&>(*this);
295  }
296 
299  const Implementation& asImpl() const
300  {
301  return static_cast<const Implementation&>(*this);
302  }
303 
304  // return true if wells are available in the reservoir
305  bool wellsActive() const { return wells_active_; }
306  // return true if wells are available on this process
307  bool localWellsActive() const { return wells_ ? (wells_->number_of_wells > 0 ) : false; }
308  // return wells object
309  const Wells& wells () const { assert( bool(wells_ != 0) ); return *wells_; }
310 
311  void
312  makeConstantState(SolutionState& state) const;
313 
314  SolutionState
315  variableState(const ReservoirState& x,
316  const WellState& xw) const;
317 
318  std::vector<V>
319  variableStateInitials(const ReservoirState& x,
320  const WellState& xw) const;
321  void
322  variableReservoirStateInitials(const ReservoirState& x,
323  std::vector<V>& vars0) const;
324  void
325  variableWellStateInitials(const WellState& xw,
326  std::vector<V>& vars0) const;
327 
328  std::vector<int>
329  variableStateIndices() const;
330 
331  std::vector<int>
332  variableWellStateIndices() const;
333 
334  SolutionState
335  variableStateExtractVars(const ReservoirState& x,
336  const std::vector<int>& indices,
337  std::vector<ADB>& vars) const;
338 
339  void
340  variableStateExtractWellsVars(const std::vector<int>& indices,
341  std::vector<ADB>& vars,
342  SolutionState& state) const;
343 
344  void
345  computeAccum(const SolutionState& state,
346  const int aix );
347 
348  void computeWellConnectionPressures(const SolutionState& state,
349  const WellState& xw);
350 
351  void
352  assembleMassBalanceEq(const SolutionState& state);
353 
354  void
355  solveWellEq(const std::vector<ADB>& mob_perfcells,
356  const std::vector<ADB>& b_perfcells,
357  SolutionState& state,
358  WellState& well_state);
359 
360  void
361  computeWellFlux(const SolutionState& state,
362  const std::vector<ADB>& mob_perfcells,
363  const std::vector<ADB>& b_perfcells,
364  V& aliveWells,
365  std::vector<ADB>& cq_s);
366 
367  void
368  updatePerfPhaseRatesAndPressures(const std::vector<ADB>& cq_s,
369  const SolutionState& state,
370  WellState& xw);
371 
372  void
373  addWellFluxEq(const std::vector<ADB>& cq_s,
374  const SolutionState& state);
375 
376  void
377  addWellContributionToMassBalanceEq(const std::vector<ADB>& cq_s,
378  const SolutionState& state,
379  const WellState& xw);
380 
381  void
382  addWellControlEq(const SolutionState& state,
383  const WellState& xw,
384  const V& aliveWells);
385 
386  void updateWellControls(WellState& xw) const;
387 
388  void updateWellState(const V& dwells,
389  WellState& well_state);
390 
391  bool getWellConvergence(const int iteration);
392 
393  bool isVFPActive() const;
394 
395  std::vector<ADB>
396  computePressures(const ADB& po,
397  const ADB& sw,
398  const ADB& so,
399  const ADB& sg) const;
400 
401  V
402  computeGasPressure(const V& po,
403  const V& sw,
404  const V& so,
405  const V& sg) const;
406 
407  std::vector<ADB>
408  computeRelPerm(const SolutionState& state) const;
409 
410  void
411  computeMassFlux(const int actph ,
412  const V& transi,
413  const ADB& kr ,
414  const ADB& p ,
415  const SolutionState& state );
416 
417  void applyThresholdPressures(ADB& dp);
418 
419  ADB
420  fluidViscosity(const int phase,
421  const ADB& p ,
422  const ADB& temp ,
423  const ADB& rs ,
424  const ADB& rv ,
425  const std::vector<PhasePresence>& cond) const;
426 
427  ADB
428  fluidReciprocFVF(const int phase,
429  const ADB& p ,
430  const ADB& temp ,
431  const ADB& rs ,
432  const ADB& rv ,
433  const std::vector<PhasePresence>& cond) const;
434 
435  ADB
436  fluidDensity(const int phase,
437  const ADB& b,
438  const ADB& rs,
439  const ADB& rv) const;
440 
441  V
442  fluidRsSat(const V& p,
443  const V& so,
444  const std::vector<int>& cells) const;
445 
446  ADB
447  fluidRsSat(const ADB& p,
448  const ADB& so,
449  const std::vector<int>& cells) const;
450 
451  V
452  fluidRvSat(const V& p,
453  const V& so,
454  const std::vector<int>& cells) const;
455 
456  ADB
457  fluidRvSat(const ADB& p,
458  const ADB& so,
459  const std::vector<int>& cells) const;
460 
461  ADB
462  poroMult(const ADB& p) const;
463 
464  ADB
465  transMult(const ADB& p) const;
466 
467  const std::vector<PhasePresence>
468  phaseCondition() const {return phaseCondition_;}
469 
470  void
471  classifyCondition(const ReservoirState& state);
472 
473 
476  void
477  updatePrimalVariableFromState(const ReservoirState& state);
478 
481  void
483 
505  double
506  convergenceReduction(const Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic>& B,
507  const Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic>& tempV,
508  const Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic>& R,
509  std::vector<double>& R_sum,
510  std::vector<double>& maxCoeff,
511  std::vector<double>& B_avg,
512  std::vector<double>& maxNormWell,
513  int nc,
514  int nw) const;
515 
516  double dpMaxRel() const { return param_.dp_max_rel_; }
517  double dsMax() const { return param_.ds_max_; }
518  double drMaxRel() const { return param_.dr_max_rel_; }
519  double maxResidualAllowed() const { return param_.max_residual_allowed_; }
520 
521  };
522 } // namespace Opm
523 
525 
526 #endif // OPM_BLACKOILMODELBASE_HEADER_INCLUDED
Definition: BlackoilModelBase.hpp:232
Definition: GeoProps.hpp:53
std::vector< double > computeResidualNorms() const
Compute the residual norms of the mass balance for each phase, the well flux, and the well equation...
Definition: BlackoilModelBase_impl.hpp:2292
void makeConstantState(SolutionState &state) const
Definition: BlackoilModelBase_impl.hpp:388
double drMaxRel() const
Definition: BlackoilModelBase.hpp:518
SolutionState variableStateExtractVars(const ReservoirState &x, const std::vector< int > &indices, std::vector< ADB > &vars) const
Definition: BlackoilModelBase_impl.hpp:566
ADB bhp
Definition: BlackoilModelBase.hpp:71
LinearisedBlackoilResidual residual_
Definition: BlackoilModelBase.hpp:279
ADB mflux
Definition: BlackoilModelBase.hpp:235
Definition: BlackoilPropsAdInterface.hpp:38
std::vector< int > primalVariable_
Definition: BlackoilModelBase.hpp:284
const DerivedGeology & geo_
Definition: BlackoilModelBase.hpp:251
Definition: LinearisedBlackoilResidual.hpp:47
ADB pressure
Definition: BlackoilModelBase.hpp:65
Eigen::Array< Scalar, Eigen::Dynamic, 1 > V
Underlying type for values.
Definition: AutoDiffBlock.hpp:98
V fluidRvSat(const V &p, const V &so, const std::vector< int > &cells) const
Definition: BlackoilModelBase_impl.hpp:2722
ADB dh
Definition: BlackoilModelBase.hpp:237
ADB qs
Definition: BlackoilModelBase.hpp:70
ADB b
Definition: BlackoilModelBase.hpp:236
Definition: AdditionalObjectDeleter.hpp:22
V isRs_
Definition: BlackoilModelBase.hpp:273
const Wells * wells_
Definition: BlackoilModelBase.hpp:253
bool getConvergence(const double dt, const int iteration)
Definition: BlackoilModelBase_impl.hpp:2420
void afterStep(const double dt, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilModelBase_impl.hpp:230
void updateState(const V &dx, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilModelBase_impl.hpp:1816
ModelTraits< Implementation >::ModelParameters ModelParameters
Definition: BlackoilModelBase.hpp:108
DefaultBlackoilSolutionState(const int np)
Definition: BlackoilModelBase.hpp:54
AutoDiffBlock< double > ADB
Definition: BlackoilModelBase.hpp:53
const Wells & wells() const
Definition: BlackoilModelBase.hpp:309
const std::vector< bool > active_
Definition: BlackoilModelBase.hpp:257
double convergenceReduction(const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &B, const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &tempV, const Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic > &R, std::vector< double > &R_sum, std::vector< double > &maxCoeff, std::vector< double > &B_avg, std::vector< double > &maxNormWell, int nc, int nw) const
Compute the reduction within the convergence check.
Definition: BlackoilModelBase_impl.hpp:2335
AutoDiffBlock< double > ADB
Definition: BlackoilModelBase.hpp:102
std::vector< std::string > material_name_
Definition: BlackoilModelBase.hpp:286
void variableStateExtractWellsVars(const std::vector< int > &indices, std::vector< ADB > &vars, SolutionState &state) const
Definition: BlackoilModelBase_impl.hpp:638
const RockCompressibility * rock_comp_props_
Definition: BlackoilModelBase.hpp:252
void computeWellConnectionPressures(const SolutionState &state, const WellState &xw)
Definition: BlackoilModelBase_impl.hpp:714
ModelTraits< Implementation >::ReservoirState ReservoirState
Definition: BlackoilModelBase.hpp:106
void updateWellState(const V &dwells, WellState &well_state)
Definition: BlackoilModelBase_impl.hpp:2040
const std::vector< int > canph_
Definition: BlackoilModelBase.hpp:259
bool wellsActive() const
Definition: BlackoilModelBase.hpp:305
Implementation & asImpl()
Definition: BlackoilModelBase.hpp:292
ADB transMult(const ADB &p) const
Definition: BlackoilModelBase_impl.hpp:2778
double maxResidualAllowed() const
Definition: BlackoilModelBase.hpp:519
std::vector< ADB > saturation
Definition: BlackoilModelBase.hpp:67
int numMaterials() const
Definition: BlackoilModelBase_impl.hpp:292
Definition: BlackoilModelBase.hpp:98
bool getWellConvergence(const int iteration)
Definition: BlackoilModelBase_impl.hpp:2542
void computeMassFlux(const int actph, const V &transi, const ADB &kr, const ADB &p, const SolutionState &state)
Definition: BlackoilModelBase_impl.hpp:2225
V pvdt_
Definition: BlackoilModelBase.hpp:285
bool localWellsActive() const
Definition: BlackoilModelBase.hpp:307
V well_perforation_pressure_diffs_
Definition: BlackoilModelBase.hpp:277
void classifyCondition(const ReservoirState &state)
Definition: BlackoilModelBase_impl.hpp:2808
const WellOps wops_
Definition: BlackoilModelBase.hpp:262
Definition: BlackoilModelBase.hpp:241
ADB::M M
Definition: BlackoilModelBase.hpp:104
const BlackoilPropsAdInterface & fluid_
Definition: BlackoilModelBase.hpp:250
ModelTraits< Implementation >::SolutionState SolutionState
Definition: BlackoilModelBase.hpp:109
V threshold_pressures_by_interior_face_
Definition: BlackoilModelBase.hpp:269
ModelParameters param_
Definition: BlackoilModelBase.hpp:266
void setThresholdPressures(const std::vector< double > &threshold_pressures_by_face)
Set threshold pressures that prevent or reduce flow. This prevents flow across faces if the potential...
Definition: BlackoilModelBase_impl.hpp:317
Struct for containing iteration variables.
Definition: BlackoilModelBase.hpp:51
std::vector< ADB > accum
Definition: BlackoilModelBase.hpp:234
Eigen::SparseMatrix< double > p2w
Definition: BlackoilModelBase.hpp:244
void updatePhaseCondFromPrimalVariable()
Update the phaseCondition_ member based on the primalVariable_ member.
Definition: BlackoilModelBase_impl.hpp:2898
std::vector< ADB > canonical_phase_pressures
Definition: BlackoilModelBase.hpp:73
std::vector< int > variableStateIndices() const
Definition: BlackoilModelBase_impl.hpp:525
BlackoilModelBase(const ModelParameters &param, const Grid &grid, const BlackoilPropsAdInterface &fluid, const DerivedGeology &geo, const RockCompressibility *rock_comp_props, const Wells *wells, const NewtonIterationBlackoilInterface &linsolver, Opm::EclipseStateConstPtr eclState, const bool has_disgas, const bool has_vapoil, const bool terminal_output)
Definition: BlackoilModelBase_impl.hpp:145
bool wells_active_
Definition: BlackoilModelBase.hpp:268
const std::vector< int > cells_
Definition: BlackoilModelBase.hpp:260
void solveWellEq(const std::vector< ADB > &mob_perfcells, const std::vector< ADB > &b_perfcells, SolutionState &state, WellState &well_state)
Definition: BlackoilModelBase_impl.hpp:1486
SolutionState variableState(const ReservoirState &x, const WellState &xw) const
Definition: BlackoilModelBase_impl.hpp:417
const std::vector< PhasePresence > phaseCondition() const
Definition: BlackoilModelBase.hpp:468
ADB rs
Definition: BlackoilModelBase.hpp:68
VFPProperties vfp_properties_
Definition: BlackoilModelBase.hpp:254
Definition: AutoDiffMatrix.hpp:43
void variableReservoirStateInitials(const ReservoirState &x, std::vector< V > &vars0) const
Definition: BlackoilModelBase_impl.hpp:453
bool use_threshold_pressure_
Definition: BlackoilModelBase.hpp:267
int linearIterationsLastSolve() const
Number of linear iterations used in last call to solveJacobianSystem().
Definition: BlackoilModelBase_impl.hpp:256
const Grid & grid_
Definition: BlackoilModelBase.hpp:249
int numPhases() const
The number of active fluid phases in the model.
Definition: BlackoilModelBase_impl.hpp:280
Interface class for (linear) solvers for the fully implicit black-oil system.
Definition: NewtonIterationBlackoilInterface.hpp:31
V well_perforation_densities_
Definition: BlackoilModelBase.hpp:276
Definition: AutoDiffHelpers.hpp:40
V isRv_
Definition: BlackoilModelBase.hpp:274
const bool has_vapoil_
Definition: BlackoilModelBase.hpp:264
ADB fluidDensity(const int phase, const ADB &b, const ADB &rs, const ADB &rv) const
Definition: BlackoilModelBase_impl.hpp:2672
std::vector< PhasePresence > phaseCondition_
Definition: BlackoilModelBase.hpp:272
HelperOps ops_
Definition: BlackoilModelBase.hpp:261
std::vector< V > variableStateInitials(const ReservoirState &x, const WellState &xw) const
Definition: BlackoilModelBase_impl.hpp:431
double dsMax() const
Definition: BlackoilModelBase.hpp:517
Eigen::SparseMatrix< double > w2p
Definition: BlackoilModelBase.hpp:243
const bool has_disgas_
Definition: BlackoilModelBase.hpp:263
ADB poroMult(const ADB &p) const
Definition: BlackoilModelBase_impl.hpp:2748
void assemble(const ReservoirState &reservoir_state, WellState &well_state, const bool initial_assembly)
Definition: BlackoilModelBase_impl.hpp:815
Eigen::Array< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > DataBlock
Definition: BlackoilModelBase.hpp:230
void computeAccum(const SolutionState &state, const int aix)
Definition: BlackoilModelBase_impl.hpp:655
V solveJacobianSystem() const
Definition: BlackoilModelBase_impl.hpp:1746
V isSg_
Definition: BlackoilModelBase.hpp:275
void updateEquationsScaling()
Update the scaling factors for mass balance equations.
Definition: BlackoilModelBase_impl.hpp:963
WellOps(const Wells *wells)
Definition: BlackoilModelBase_impl.hpp:352
std::vector< ReservoirResidualQuant > rq_
Definition: BlackoilModelBase.hpp:271
Definition: VFPProperties.hpp:37
bool terminal_output_
Whether we print something to std::cout.
Definition: BlackoilModelBase.hpp:282
int sizeNonLinear() const
The size (number of unknowns) of the nonlinear system of equations.
Definition: BlackoilModelBase_impl.hpp:244
V fluidRsSat(const V &p, const V &so, const std::vector< int > &cells) const
Definition: BlackoilModelBase_impl.hpp:2696
void prepareStep(const double dt, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilModelBase_impl.hpp:213
void updatePrimalVariableFromState(const ReservoirState &state)
Definition: BlackoilModelBase_impl.hpp:2848
double dpMaxRel() const
Definition: BlackoilModelBase.hpp:516
ADB fluidReciprocFVF(const int phase, const ADB &p, const ADB &temp, const ADB &rs, const ADB &rv, const std::vector< PhasePresence > &cond) const
Definition: BlackoilModelBase_impl.hpp:2647
void addWellFluxEq(const std::vector< ADB > &cq_s, const SolutionState &state)
Definition: BlackoilModelBase_impl.hpp:1177
const Implementation & asImpl() const
Definition: BlackoilModelBase.hpp:299
void updateWellControls(WellState &xw) const
Definition: BlackoilModelBase_impl.hpp:1360
std::vector< ADB > computeRelPerm(const SolutionState &state) const
Definition: BlackoilModelBase_impl.hpp:2140
void computeWellFlux(const SolutionState &state, const std::vector< ADB > &mob_perfcells, const std::vector< ADB > &b_perfcells, V &aliveWells, std::vector< ADB > &cq_s)
Definition: BlackoilModelBase_impl.hpp:1003
void assembleMassBalanceEq(const SolutionState &state)
Definition: BlackoilModelBase_impl.hpp:899
std::vector< int > variableWellStateIndices() const
Definition: BlackoilModelBase_impl.hpp:548
ADB rv
Definition: BlackoilModelBase.hpp:69
ADB mob
Definition: BlackoilModelBase.hpp:238
V computeGasPressure(const V &po, const V &sw, const V &so, const V &sg) const
Definition: BlackoilModelBase_impl.hpp:2206
bool terminalOutputEnabled() const
Return true if output to cout is wanted.
Definition: BlackoilModelBase_impl.hpp:268
void addWellContributionToMassBalanceEq(const std::vector< ADB > &cq_s, const SolutionState &state, const WellState &xw)
Definition: BlackoilModelBase_impl.hpp:983
ADB fluidViscosity(const int phase, const ADB &p, const ADB &temp, const ADB &rs, const ADB &rv, const std::vector< PhasePresence > &cond) const
Definition: BlackoilModelBase_impl.hpp:2622
ReservoirResidualQuant()
Definition: BlackoilModelBase_impl.hpp:337
const std::string & materialName(int material_index) const
Definition: BlackoilModelBase_impl.hpp:304
ADB::V V
Definition: BlackoilModelBase.hpp:103
void updatePerfPhaseRatesAndPressures(const std::vector< ADB > &cq_s, const SolutionState &state, WellState &xw)
Definition: BlackoilModelBase_impl.hpp:1152
bool isVFPActive() const
Definition: BlackoilModelBase_impl.hpp:1328
Definition: BlackoilModelBase.hpp:83
void applyThresholdPressures(ADB &dp)
Definition: BlackoilModelBase_impl.hpp:2260
void addWellControlEq(const SolutionState &state, const WellState &xw, const V &aliveWells)
Definition: BlackoilModelBase_impl.hpp:1580
ModelTraits< Implementation >::WellState WellState
Definition: BlackoilModelBase.hpp:107
ADB temperature
Definition: BlackoilModelBase.hpp:66
void variableWellStateInitials(const WellState &xw, std::vector< V > &vars0) const
Definition: BlackoilModelBase_impl.hpp:491
const NewtonIterationBlackoilInterface & linsolver_
Definition: BlackoilModelBase.hpp:255
std::vector< ADB > computePressures(const ADB &po, const ADB &sw, const ADB &so, const ADB &sg) const
Definition: BlackoilModelBase_impl.hpp:2170