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::WellGroupHelperType;
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 WellGroupHelperType& wgHelper,
95 WellStateType& well_state,
96 DeferredLogger& deferred_logger) const override;
97
99 void scaleSegmentRatesAndPressure(WellStateType& well_state) const override;
100
103 const WellStateType& well_state,
104 const std::vector<Scalar>& B_avg,
105 DeferredLogger& deferred_logger,
106 const bool relax_tolerance) const override;
107
109 void apply(const BVector& x, BVector& Ax) const override;
111 void apply(BVector& r) const override;
112
116 const BVector& x,
117 WellStateType& well_state,
118 DeferredLogger& deferred_logger) override;
119
121 void computeWellPotentials(const Simulator& simulator,
122 const WellStateType& well_state,
123 const WellGroupHelperType& wgHelper,
124 std::vector<Scalar>& well_potentials,
125 DeferredLogger& deferred_logger) override;
126
127 void updatePrimaryVariables(const Simulator& simulator,
128 const WellStateType& well_state,
129 DeferredLogger& deferred_logger) override;
130
131 void solveEqAndUpdateWellState(const Simulator& simulator,
132 WellStateType& well_state,
133 DeferredLogger& deferred_logger) override; // const?
134
135 void calculateExplicitQuantities(const Simulator& simulator,
136 const WellStateType& well_state,
137 DeferredLogger& deferred_logger) override; // should be const?
138
139 void updateIPRImplicit(const Simulator& simulator,
140 const WellGroupHelperType& wgHelper,
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 WellGroupHelperType& wgHelper,
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 WellStateType& well_state,
189 DeferredLogger& deferred_logger,
190 const Scalar relaxation_factor = 1.0);
191
192 // computing the accumulation term for later use in well mass equations
193 void computeInitialSegmentFluids(const Simulator& simulator, DeferredLogger& deferred_logger);
194
195 // compute the pressure difference between the perforation and cell center
196 void computePerfCellPressDiffs(const Simulator& simulator);
197
198 template<class Value>
199 void computePerfRate(const IntensiveQuantities& int_quants,
200 const std::vector<Value>& mob_perfcells,
201 const std::vector<Value>& Tw,
202 const int seg,
203 const int perf,
204 const Value& segment_pressure,
205 const bool& allow_cf,
206 std::vector<Value>& cq_s,
207 Value& perf_press,
208 PerforationRates<Scalar>& perf_rates,
209 DeferredLogger& deferred_logger) const;
210
211 template<class Value>
212 void computePerfRate(const Value& pressure_cell,
213 const Value& rs,
214 const Value& rv,
215 const std::vector<Value>& b_perfcells,
216 const std::vector<Value>& mob_perfcells,
217 const std::vector<Value>& Tw,
218 const int perf,
219 const Value& segment_pressure,
220 const Value& segment_density,
221 const bool& allow_cf,
222 const std::vector<Value>& cmix_s,
223 std::vector<Value>& cq_s,
224 Value& perf_press,
225 PerforationRates<Scalar>& perf_rates,
226 DeferredLogger& deferred_logger) const;
227
228 // compute the fluid properties, such as densities, viscosities, and so on, in the segments
229 // They will be treated implicitly, so they need to be of Evaluation type
230 void computeSegmentFluidProperties(const Simulator& simulator,
231 DeferredLogger& deferred_logger);
232
233 // get the transmissibility multiplier for specific perforation
234 template<class Value>
235 void getTransMult(Value& trans_mult,
236 const Simulator& simulator,
237 const int cell_indx) const;
238
239 // get the mobility for specific perforation
240 template<class Value>
241 void getMobility(const Simulator& simulator,
242 const int local_perf_index,
243 std::vector<Value>& mob,
244 DeferredLogger& deferred_logger) const;
245
246 void computeWellRatesAtBhpLimit(const Simulator& simulator,
247 const WellGroupHelperType& wgHelper,
248 std::vector<Scalar>& well_flux,
249 DeferredLogger& deferred_logger) const;
250
251 void computeWellRatesWithBhp(const Simulator& simulator,
252 const Scalar& bhp,
253 std::vector<Scalar>& well_flux,
254 DeferredLogger& deferred_logger) const override;
255
256 void computeWellRatesWithBhpIterations(const Simulator& simulator,
257 const Scalar& bhp,
258 const WellGroupHelperType& wgHelper,
259 std::vector<Scalar>& well_flux,
260 DeferredLogger& deferred_logger) const override;
261
262 std::vector<Scalar>
264 const Simulator& simulator,
265 const WellGroupHelperType& wgHelper,
266 DeferredLogger& deferred_logger) const;
267
268 bool computeWellPotentialsImplicit(const Simulator& simulator,
269 const WellGroupHelperType& wgHelper,
270 std::vector<Scalar>& well_potentials,
271 DeferredLogger& deferred_logger) const;
272
273 Scalar getRefDensity() const override;
274
275 bool iterateWellEqWithControl(const Simulator& simulator,
276 const double dt,
277 const Well::InjectionControls& inj_controls,
278 const Well::ProductionControls& prod_controls,
279 const WellGroupHelperType& wgHelper,
280 WellStateType& well_state,
281 DeferredLogger& deferred_logger) override;
282
283 bool iterateWellEqWithSwitching(const Simulator& simulator,
284 const double dt,
285 const Well::InjectionControls& inj_controls,
286 const Well::ProductionControls& prod_controls,
287 const WellGroupHelperType& wgHelper,
288 WellStateType& well_state,
289 DeferredLogger& deferred_logger,
290 const bool fixed_control,
291 const bool fixed_status,
292 const bool solving_with_zero_rate) override;
293
294 void assembleWellEqWithoutIteration(const Simulator& simulator,
295 const WellGroupHelperType& wgHelper,
296 const double dt,
297 const Well::InjectionControls& inj_controls,
298 const Well::ProductionControls& prod_controls,
299 WellStateType& well_state,
300 DeferredLogger& deferred_logger,
301 const bool solving_with_zero_rate) override;
302
303 void updateWaterThroughput(const double dt, WellStateType& well_state) const override;
304
306 const int seg_idx,
307 DeferredLogger& deferred_logger) const;
308
309 // turn on crossflow to avoid singular well equations
310 // when the well is banned from cross-flow and the BHP is not properly initialized,
311 // we turn on crossflow to avoid singular well equations. It can result in wrong-signed
312 // well rates, it can cause problem for THP calculation
313 // TODO: looking for better alternative to avoid wrong-signed well rates
314 bool openCrossFlowAvoidSingularity(const Simulator& simulator) const;
315
316 // for a well, when all drawdown are in the wrong direction, then this well will not
317 // be able to produce/inject .
318 bool allDrawDownWrongDirection(const Simulator& simulator) const;
319
320 std::optional<Scalar>
321 computeBhpAtThpLimitProd(const WellStateType& well_state,
322 const Simulator& ebos_simulator,
323 const WellGroupHelperType& wgHelper,
324 const SummaryState& summary_state,
325 DeferredLogger& deferred_logger) const;
326
327 std::optional<Scalar>
328 computeBhpAtThpLimitInj(const Simulator& ebos_simulator,
329 const WellGroupHelperType& wgHelper,
330 const SummaryState& summary_state,
331 DeferredLogger& deferred_logger) const;
332
333 Scalar maxPerfPress(const Simulator& simulator) const;
334
335 // check whether the well is operable under BHP limit with current reservoir condition
336 void checkOperabilityUnderBHPLimit(const WellStateType& well_state,
337 const Simulator& ebos_simulator,
338 DeferredLogger& deferred_logger) override;
339
340 // check whether the well is operable under THP limit with current reservoir condition
341 void checkOperabilityUnderTHPLimit(const Simulator& ebos_simulator,
342 const WellStateType& well_state,
343 const WellGroupHelperType& wgHelper,
344 DeferredLogger& deferred_logger) override;
345
346 // updating the inflow based on the current reservoir condition
347 void updateIPR(const Simulator& ebos_simulator,
348 DeferredLogger& deferred_logger) const override;
349
350 FSInfo getFirstPerforationFluidStateInfo(const Simulator& simulator) const;
351 };
352
353}
354
356
357#endif // OPM_MULTISEGMENTWELL_HEADER_INCLUDED
Definition: ConvergenceReport.hpp:38
Definition: DeferredLogger.hpp:57
Definition: MultisegmentWellEval.hpp:48
MultisegmentWellEquations< Scalar, IndexTraits, numWellEq, Indices::numEq > Equations
Definition: MultisegmentWellEval.hpp:57
Definition: MultisegmentWell.hpp:38
std::optional< Scalar > computeBhpAtThpLimitProd(const WellStateType &well_state, const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2117
void updateWellState(const Simulator &simulator, const BVectorWell &dwells, WellStateType &well_state, DeferredLogger &deferred_logger, const Scalar relaxation_factor=1.0)
Definition: MultisegmentWell_impl.hpp:722
void updateWaterThroughput(const double dt, WellStateType &well_state) const override
Definition: MultisegmentWell_impl.hpp:2084
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:894
Scalar connectionDensity(const int globalConnIdx, const int openConnIdx) const override
Definition: MultisegmentWell_impl.hpp:859
void addWellContributions(SparseMatrixAdapter &jacobian) const override
Definition: MultisegmentWell_impl.hpp:881
bool iterateWellEqWithControl(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1525
void getTransMult(Value &trans_mult, const Simulator &simulator, const int cell_indx) const
Definition: MultisegmentWell_impl.hpp:1167
void solveEqAndUpdateWellState(const Simulator &simulator, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:617
Scalar getRefDensity() const override
Definition: MultisegmentWell_impl.hpp:1227
std::optional< Scalar > computeBhpAtThpLimitProdWithAlq(const Simulator &simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, const Scalar alq_value, DeferredLogger &deferred_logger, bool iterate_if_no_solution) const override
Definition: MultisegmentWell_impl.hpp:2137
std::vector< Scalar > getPrimaryVars() const override
Definition: MultisegmentWell_impl.hpp:2323
void updateIPRImplicit(const Simulator &simulator, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1402
void assembleWellEqWithoutIteration(const Simulator &simulator, const WellGroupHelperType &wgHelper, 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:1829
void computeWellRatesWithBhpIterations(const Simulator &simulator, const Scalar &bhp, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:410
bool iterateWellEqWithSwitching(const Simulator &simulator, const double dt, const Well::InjectionControls &inj_controls, const Well::ProductionControls &prod_controls, const WellGroupHelperType &wgHelper, 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:1650
void computeWellPotentials(const Simulator &simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) override
computing the well potentials for group control
Definition: MultisegmentWell_impl.hpp:297
std::vector< Scalar > computeCurrentWellRates(const Simulator &simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:2277
void apply(const BVector &x, BVector &Ax) const override
Ax = Ax - C D^-1 B x.
Definition: MultisegmentWell_impl.hpp:230
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:1235
void getMobility(const Simulator &simulator, const int local_perf_index, std::vector< Value > &mob, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:1187
bool openCrossFlowAvoidSingularity(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2022
void computeSegmentFluidProperties(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:1136
bool computeWellPotentialsImplicit(const Simulator &simulator, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_potentials, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:539
int setPrimaryVars(typename std::vector< Scalar >::const_iterator it) override
Definition: MultisegmentWell_impl.hpp:2341
void computeWellRatesWithBhp(const Simulator &simulator, const Scalar &bhp, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:360
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:2031
int debug_cost_counter_
Definition: MultisegmentWell.hpp:183
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: MultisegmentWell_impl.hpp:204
void updateProductivityIndex(const Simulator &simulator, const WellProdIndexCalculator< Scalar > &wellPICalc, WellStateType &well_state, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:777
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:1063
std::optional< Scalar > computeBhpAtThpLimitInj(const Simulator &ebos_simulator, const WellGroupHelperType &wgHelper, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2194
void computePerfCellPressDiffs(const Simulator &simulator)
Definition: MultisegmentWell_impl.hpp:646
void computeWellRatesAtBhpLimit(const Simulator &simulator, const WellGroupHelperType &wgHelper, std::vector< Scalar > &well_flux, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:343
void recoverWellSolutionAndUpdateWellState(const Simulator &simulator, const BVector &x, WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:265
std::vector< Scalar > computeWellPotentialWithTHP(const WellStateType &well_state, const Simulator &simulator, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:487
void updateIPR(const Simulator &ebos_simulator, DeferredLogger &deferred_logger) const override
Definition: MultisegmentWell_impl.hpp:1301
FSInfo getFirstPerforationFluidStateInfo(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2357
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
Scalar maxPerfPress(const Simulator &simulator) const
Definition: MultisegmentWell_impl.hpp:2249
void computeInitialSegmentFluids(const Simulator &simulator, DeferredLogger &deferred_logger)
Definition: MultisegmentWell_impl.hpp:703
std::vector< std::vector< Scalar > > segment_fluid_initial_
Definition: MultisegmentWell.hpp:181
void checkOperabilityUnderTHPLimit(const Simulator &ebos_simulator, const WellStateType &well_state, const WellGroupHelperType &wgHelper, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:1472
void updatePrimaryVariables(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:157
EvalWell getSegmentSurfaceVolume(const Simulator &simulator, const int seg_idx, DeferredLogger &deferred_logger) const
Definition: MultisegmentWell_impl.hpp:2095
void calculateExplicitQuantities(const Simulator &simulator, const WellStateType &well_state, DeferredLogger &deferred_logger) override
Definition: MultisegmentWell_impl.hpp:761
void updateWellStateWithTarget(const Simulator &simulator, const WellGroupHelperType &wgHelper, WellStateType &well_state, DeferredLogger &deferred_logger) const override
updating the well state based the current control mode
Definition: MultisegmentWell_impl.hpp:184
Class encapsulating some information about parallel wells.
Definition: ParallelWellInfo.hpp:198
Definition: WellGroupHelper.hpp:59
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
WellGroupHelper< Scalar, IndexTraits > WellGroupHelperType
Definition: WellInterface.hpp:102
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
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