MultisegmentWell.hpp
Go to the documentation of this file.
1/*
2 Copyright 2017 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2017 Statoil ASA.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21
22#ifndef OPM_MULTISEGMENTWELL_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_HEADER_INCLUDED
24
26
29
30namespace Opm {
31
32 class DeferredLogger;
33
34 template<typename TypeTag>
35 class MultisegmentWell : public WellInterface<TypeTag>
36 , public MultisegmentWellEval<GetPropType<TypeTag, Properties::FluidSystem>,
37 GetPropType<TypeTag, Properties::Indices>>
38 {
39 public:
43
44 using typename Base::Simulator;
45 using typename Base::IntensiveQuantities;
46 using typename Base::FluidSystem;
47 using typename Base::IndexTraits;
48 using typename Base::ModelParameters;
49 using typename Base::MaterialLaw;
50 using typename Base::Indices;
51 using typename Base::RateConverterType;
52 using typename Base::SparseMatrixAdapter;
53 using typename Base::FluidState;
54 using typename Base::WellStateType;
55 using typename Base::GroupStateHelperType;
56
59 using Base::Water;
60 using Base::Oil;
61 using Base::Gas;
62
63 using typename Base::Scalar;
64
66 using typename Base::BVector;
67 using typename Base::Eval;
68
69 using typename MSWEval::Equations;
70 using typename MSWEval::EvalWell;
71 using typename MSWEval::BVectorWell;
72 using MSWEval::SPres;
73 using typename Base::PressureMatrix;
74 using FSInfo = std::tuple<Scalar,typename std::decay<decltype(std::declval<decltype(std::declval<const Simulator&>().model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type>;
75
76 MultisegmentWell(const Well& well,
77 const ParallelWellInfo<Scalar>& pw_info,
78 const int time_step,
79 const ModelParameters& param,
80 const RateConverterType& rate_converter,
81 const int pvtRegionIdx,
82 const int num_conservation_quantities,
83 const int num_phases,
84 const int index_of_well,
85 const std::vector<PerforationData<Scalar>>& perf_data);
86
87 void init(const std::vector<Scalar>& depth_arg,
88 const Scalar gravity_arg,
89 const std::vector<Scalar>& B_avg,
90 const bool changed_to_open_this_step) override;
91
93 void updateWellStateWithTarget(const Simulator& simulator,
94 const GroupStateHelperType& groupStateHelper,
95 WellStateType& well_state,
96 DeferredLogger& deferred_logger) const override;
97
99 void scaleSegmentRatesAndPressure(WellStateType& well_state) const override;
100
103 const std::vector<Scalar>& B_avg,
104 DeferredLogger& deferred_logger,
105 const bool relax_tolerance) const override;
106
108 void apply(const BVector& x, BVector& Ax) const override;
110 void apply(BVector& r) const override;
111
115 const BVector& x,
116 const GroupStateHelperType& groupStateHelper,
117 WellStateType& well_state,
118 DeferredLogger& deferred_logger) override;
119
121 void computeWellPotentials(const Simulator& simulator,
122 const WellStateType& well_state,
123 const GroupStateHelperType& groupStateHelper,
124 std::vector<Scalar>& well_potentials,
125 DeferredLogger& deferred_logger) override;
126
127 void updatePrimaryVariables(const GroupStateHelperType& groupStateHelper,
128 DeferredLogger& deferred_logger) override;
129
130 void solveEqAndUpdateWellState(const Simulator& simulator,
131 const GroupStateHelperType& groupStateHelper,
132 WellStateType& well_state,
133 DeferredLogger& deferred_logger) override; // const?
134
135 void calculateExplicitQuantities(const Simulator& simulator,
136 const GroupStateHelperType& groupStateHelper,
137 DeferredLogger& deferred_logger) override; // should be const?
138
139 void updateIPRImplicit(const Simulator& simulator,
140 const GroupStateHelperType& groupStateHelper,
141 WellStateType& well_state,
142 DeferredLogger& deferred_logger) override;
143
144 void updateProductivityIndex(const Simulator& simulator,
145 const WellProdIndexCalculator<Scalar>& wellPICalc,
146 WellStateType& well_state,
147 DeferredLogger& deferred_logger) const override;
148
149 Scalar connectionDensity(const int globalConnIdx,
150 const int openConnIdx) const override;
151
152 void addWellContributions(SparseMatrixAdapter& jacobian) const override;
153
155 const BVector& x,
156 const int pressureVarIndex,
157 const bool use_well_weights,
158 const WellStateType& well_state) const override;
159
160 std::vector<Scalar>
161 computeCurrentWellRates(const Simulator& simulator,
162 DeferredLogger& deferred_logger) const override;
163
164 std::optional<Scalar>
166 const GroupStateHelperType& groupStateHelper,
167 const SummaryState& summary_state,
168 const Scalar alq_value,
169 DeferredLogger& deferred_logger,
170 bool iterate_if_no_solution) const override;
171
172 std::vector<Scalar> getPrimaryVars() const override;
173
174 int setPrimaryVars(typename std::vector<Scalar>::const_iterator it) override;
175
176 protected:
177 // regularize msw equation
179
180 // the intial amount of fluids in each segment under surface condition
181 std::vector<std::vector<Scalar> > segment_fluid_initial_;
182
183 mutable int debug_cost_counter_ = 0;
184
185 // updating the well_state based on well solution dwells
186 void updateWellState(const Simulator& simulator,
187 const BVectorWell& dwells,
188 const GroupStateHelperType& groupStateHelper,
189 WellStateType& well_state,
190 DeferredLogger& deferred_logger,
191 const Scalar relaxation_factor = 1.0);
192
193 // computing the accumulation term for later use in well mass equations
194 void computeInitialSegmentFluids(const Simulator& simulator, DeferredLogger& deferred_logger);
195
196 // compute the pressure difference between the perforation and cell center
197 void computePerfCellPressDiffs(const Simulator& simulator);
198
199 template<class Value>
200 void computePerfRate(const IntensiveQuantities& int_quants,
201 const std::vector<Value>& mob_perfcells,
202 const std::vector<Value>& Tw,
203 const int seg,
204 const int perf,
205 const Value& segment_pressure,
206 const bool& allow_cf,
207 std::vector<Value>& cq_s,
208 Value& perf_press,
209 PerforationRates<Scalar>& perf_rates,
210 DeferredLogger& deferred_logger) const;
211
212 template<class Value>
213 void computePerfRate(const Value& pressure_cell,
214 const Value& rs,
215 const Value& rv,
216 const std::vector<Value>& b_perfcells,
217 const std::vector<Value>& mob_perfcells,
218 const std::vector<Value>& Tw,
219 const int perf,
220 const Value& segment_pressure,
221 const Value& segment_density,
222 const bool& allow_cf,
223 const std::vector<Value>& cmix_s,
224 std::vector<Value>& cq_s,
225 Value& perf_press,
226 PerforationRates<Scalar>& perf_rates,
227 DeferredLogger& deferred_logger) const;
228
229 // compute the fluid properties, such as densities, viscosities, and so on, in the segments
230 // They will be treated implicitly, so they need to be of Evaluation type
231 void computeSegmentFluidProperties(const Simulator& simulator,
232 DeferredLogger& deferred_logger);
233
234 // get the transmissibility multiplier for specific perforation
235 template<class Value>
236 void getTransMult(Value& trans_mult,
237 const Simulator& simulator,
238 const int cell_indx) const;
239
240 // get the mobility for specific perforation
241 template<class Value>
242 void getMobility(const Simulator& simulator,
243 const int local_perf_index,
244 std::vector<Value>& mob,
245 DeferredLogger& deferred_logger) const;
246
247 void computeWellRatesAtBhpLimit(const Simulator& simulator,
248 const GroupStateHelperType& groupStateHelper,
249 std::vector<Scalar>& well_flux,
250 DeferredLogger& deferred_logger) const;
251
252 void computeWellRatesWithBhp(const Simulator& simulator,
253 const Scalar& bhp,
254 std::vector<Scalar>& well_flux,
255 DeferredLogger& deferred_logger) const override;
256
257 void computeWellRatesWithBhpIterations(const Simulator& simulator,
258 const Scalar& bhp,
259 const GroupStateHelperType& groupStateHelper,
260 std::vector<Scalar>& well_flux,
261 DeferredLogger& deferred_logger) const override;
262
263 std::vector<Scalar>
265 const Simulator& simulator,
266 const GroupStateHelperType& groupStateHelper,
267 DeferredLogger& deferred_logger) const;
268
269 bool computeWellPotentialsImplicit(const Simulator& simulator,
270 const GroupStateHelperType& groupStateHelper,
271 std::vector<Scalar>& well_potentials,
272 DeferredLogger& deferred_logger) const;
273
274 Scalar getRefDensity() const override;
275
276 bool iterateWellEqWithControl(const Simulator& simulator,
277 const double dt,
278 const Well::InjectionControls& inj_controls,
279 const Well::ProductionControls& prod_controls,
280 const GroupStateHelperType& groupStateHelper,
281 WellStateType& well_state,
282 DeferredLogger& deferred_logger) override;
283
284 bool iterateWellEqWithSwitching(const Simulator& simulator,
285 const double dt,
286 const Well::InjectionControls& inj_controls,
287 const Well::ProductionControls& prod_controls,
288 const GroupStateHelperType& groupStateHelper,
289 WellStateType& well_state,
290 DeferredLogger& deferred_logger,
291 const bool fixed_control,
292 const bool fixed_status,
293 const bool solving_with_zero_rate) override;
294
295 void assembleWellEqWithoutIteration(const Simulator& simulator,
296 const GroupStateHelperType& groupStateHelper,
297 const double dt,
298 const Well::InjectionControls& inj_controls,
299 const Well::ProductionControls& prod_controls,
300 WellStateType& well_state,
301 DeferredLogger& deferred_logger,
302 const bool solving_with_zero_rate) override;
303
304 void updateWaterThroughput(const double dt, WellStateType& well_state) const override;
305
307 const int seg_idx,
308 DeferredLogger& deferred_logger) const;
309
310 // turn on crossflow to avoid singular well equations
311 // when the well is banned from cross-flow and the BHP is not properly initialized,
312 // we turn on crossflow to avoid singular well equations. It can result in wrong-signed
313 // well rates, it can cause problem for THP calculation
314 // TODO: looking for better alternative to avoid wrong-signed well rates
315 bool openCrossFlowAvoidSingularity(const Simulator& simulator) const;
316
317 // for a well, when all drawdown are in the wrong direction, then this well will not
318 // be able to produce/inject .
319 bool allDrawDownWrongDirection(const Simulator& simulator) const;
320
321 std::optional<Scalar>
322 computeBhpAtThpLimitProd(const WellStateType& well_state,
323 const Simulator& ebos_simulator,
324 const GroupStateHelperType& groupStateHelper,
325 const SummaryState& summary_state,
326 DeferredLogger& deferred_logger) const;
327
328 std::optional<Scalar>
329 computeBhpAtThpLimitInj(const Simulator& ebos_simulator,
330 const GroupStateHelperType& groupStateHelper,
331 const SummaryState& summary_state,
332 DeferredLogger& deferred_logger) const;
333
334 Scalar maxPerfPress(const Simulator& simulator) const override;
335
336 // check whether the well is operable under BHP limit with current reservoir condition
337 void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
338 const Simulator& ebos_simulator,
339 DeferredLogger& deferred_logger) override;
340
341 // check whether the well is operable under THP limit with current reservoir condition
342 void checkOperabilityUnderTHPLimit(const Simulator& ebos_simulator,
343 const WellStateType& well_state,
344 const GroupStateHelperType& groupStateHelper,
345 DeferredLogger& deferred_logger) override;
346
347 // updating the inflow based on the current reservoir condition
348 void updateIPR(const Simulator& ebos_simulator,
349 DeferredLogger& deferred_logger) const override;
350
351 FSInfo getFirstPerforationFluidStateInfo(const Simulator& simulator) const;
352 };
353
354}
355
357
358#endif // OPM_MULTISEGMENTWELL_HEADER_INCLUDED
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: GroupStateHelper.hpp:53
Definition: MultisegmentWellEval.hpp:48
MultisegmentWellEquations< Scalar, IndexTraits, numWellEq, Indices::numEq > Equations
Definition: MultisegmentWellEval.hpp:57
Definition: MultisegmentWell.hpp:38
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: MultisegmentWell_impl.hpp:2086
void addWellPressureEquations(PressureMatrix &mat, const BVector &x, const int pressureVarIndex, const bool use_well_weights, const WellStateType &well_state) const override
Definition: MultisegmentWell_impl.hpp:897
void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:411
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:862
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:884
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_indx) const
Definition: MultisegmentWell_impl.hpp:1170
void updateWellStateWithTarget(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) const override
updating the well state based the current control mode
Definition: MultisegmentWell_impl.hpp:184
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: MultisegmentWell_impl.hpp:298
void assembleWellEqWithoutIteration(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, WellStateType &well_state, DeferredLogger &deferred_logger, const bool solving_with_zero_rate) override
Definition: MultisegmentWell_impl.hpp:1832
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1230
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2325
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2196
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const override
Definition: MultisegmentWell_impl.hpp:2139
void solveEqAndUpdateWellState(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:618
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:265
std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2279
void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: MultisegmentWell_impl.hpp:230
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:724
MultisegmentWell(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)
Definition: MultisegmentWell_impl.hpp:62
void checkOperabilityUnderBHPLimit(const WellStateType &well_state, const Simulator &ebos_simulator, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1238
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1190
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger, const bool fixed_control, const bool fixed_status, const bool solving_with_zero_rate) override
Definition: MultisegmentWell_impl.hpp:1653
ConvergenceReport getWellConvergence(const GroupStateHelperType &groupStateHelper, 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: MultisegmentWell_impl.hpp:204
std::vector< Scalar > computeWellPotentialWithTHP(const WellStateType &well_state, const Simulator &simulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:488
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2024
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1139
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: MultisegmentWell_impl.hpp:2343
void computeWellRatesWithBhp(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:361
void scaleSegmentRatesAndPressure(WellStateType &well_state) const override
updating the segment pressure and rates based the current bhp and well rates
Definition: MultisegmentWell_impl.hpp:173
bool allDrawDownWrongDirection(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2033
int debug_cost_counter_
Definition: MultisegmentWell.hpp:183
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:780
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1528
Scalar maxPerfPress(const Simulator &simulator) const override
Definition: MultisegmentWell_impl.hpp:2251
void computePerfRate(const IntensiveQuantities &int_quants, const std::vector< Value > &mob_perfcells, const std::vector< Value > &Tw, const int seg, const int perf, const Value &segment_pressure, const bool &allow_cf, std::vector< Value > &cq_s, Value &perf_press, PerforationRates< Scalar > &perf_rates, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1066
void calculateExplicitQuantities(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:764
void computeWellRatesAtBhpLimit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:344
bool computeWellPotentialsImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:540
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:648
void updateIPRImplicit(const Simulator &simulator, const GroupStateHelperType &groupStateHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1405
void updateIPR(const Simulator &ebos_simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:1304
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellStateType &well_state, const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1475
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &ebos_simulator, const GroupStateHelperType &groupStateHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2119
FSInfo getFirstPerforationFluidStateInfo(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2359
std::tuple< Scalar, typename std::decay< decltype(std::declval< decltype(std::declval< const Simulator & >().model().intensiveQuantities(0, 0).fluidState())>().saltConcentration())>::type > FSInfo
Definition: MultisegmentWell.hpp:74
void init(const std::vector< Scalar > &depth_arg, const Scalar gravity_arg, const std::vector< Scalar > &B_avg, const bool changed_to_open_this_step) override
Definition: MultisegmentWell_impl.hpp:122
bool regularize_
Definition: MultisegmentWell.hpp:178
void computeInitialSegmentFluids(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:705
std::vector< std::vector< Scalar > > segment_fluid_initial_
Definition: MultisegmentWell.hpp:181
EvalWell getSegmentSurfaceVolume(const Simulator &simulator, const int seg_idx, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2097
void updatePrimaryVariables(const GroupStateHelperType &groupStateHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:157
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
static constexpr int Oil
Definition: WellInterfaceFluidSystem.hpp:65
static constexpr int Water
Definition: WellInterfaceFluidSystem.hpp:64
static constexpr int Gas
Definition: WellInterfaceFluidSystem.hpp:66
int pvtRegionIdx() const
Definition: WellInterfaceGeneric.hpp:130
Definition: WellInterface.hpp:77
BlackOilFluidState< Eval, FluidSystem, energyModuleType==EnergyModules::ConstantTemperature,(energyModuleType==EnergyModules::FullyImplicitThermal||energyModuleType==EnergyModules::SequentialImplicitThermal), Indices::compositionSwitchIdx >=0, has_watVapor, has_brine, has_saltPrecip, has_disgas_in_water, Indices::numPhases > FluidState
Definition: WellInterface.hpp:139
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: WellInterface.hpp:82
typename WellInterfaceFluidSystem< FluidSystem >::RateConverterType RateConverterType
Definition: WellInterface.hpp:105
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > PressureMatrix
Definition: WellInterface.hpp:98
GetPropType< TypeTag, Properties::IntensiveQuantities > IntensiveQuantities
Definition: WellInterface.hpp:87
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: WellInterface.hpp:83
GetPropType< TypeTag, Properties::MaterialLaw > MaterialLaw
Definition: WellInterface.hpp:88
typename FluidSystem::IndexTraitsType IndexTraits
Definition: WellInterface.hpp:85
Dune::BlockVector< VectorBlockType > BVector
Definition: WellInterface.hpp:97
static constexpr bool has_polymer
Definition: WellInterface.hpp:115
typename Base::ModelParameters ModelParameters
Definition: WellInterface.hpp:111
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: WellInterface.hpp:84
static constexpr bool has_solvent
Definition: WellInterface.hpp:113
typename Base::Eval Eval
Definition: WellInterface.hpp:96
WellState< Scalar, IndexTraits > WellStateType
Definition: WellInterface.hpp:100
GetPropType< TypeTag, Properties::Indices > Indices
Definition: WellInterface.hpp:86
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: WellInterface.hpp:89
GroupStateHelper< Scalar, IndexTraits > GroupStateHelperType
Definition: WellInterface.hpp:102
Definition: WellProdIndexCalculator.hpp:37
Definition: WellState.hpp:66
Defines the common properties required by the porous medium multi-phase models.
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