22#ifndef OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_PRIMARY_VARIABLES_HEADER_INCLUDED
25#include <opm/material/densead/Evaluation.hpp>
27#include <opm/simulators/wells/MultisegmentWellEquations.hpp>
36template<
class Scalar>
class MultisegmentWellGeneric;
37template<
class Flu
idSystem,
class Indices,
class Scalar>
class WellInterfaceIndices;
40template<
class Flu
idSystem,
class Indices,
class Scalar>
57 static constexpr bool has_water = (Indices::waterSwitchIdx >= 0);
58 static constexpr bool has_gas = (Indices::compositionSwitchIdx >= 0);
59 static constexpr bool has_oil = (Indices::numPhases - has_gas - has_water) > 0;
64 static constexpr bool has_wfrac_variable = has_water && Indices::numPhases > 1;
65 static constexpr bool has_gfrac_variable = has_gas && has_oil;
67 static constexpr int WQTotal = 0;
68 static constexpr int WFrac = has_wfrac_variable ? 1 : -1000;
69 static constexpr int GFrac = has_gfrac_variable ? has_wfrac_variable + 1 : -1000;
70 static constexpr int SPres = has_wfrac_variable + has_gfrac_variable + 1;
73 static constexpr int numWellEq = Indices::numPhases + 1;
75 using EvalWell = DenseAd::Evaluation<double, Indices::numEq + numWellEq>;
78 using BVectorWell =
typename Equations::BVectorWell;
85 void resize(
const int numSegments);
91 void update(
const WellState& well_state,
const bool stop_or_zero_rate_target);
95 const double relaxation_factor,
97 const bool stop_or_zero_rate_target,
98 const double max_pressure_change);
103 const bool stop_or_zero_rate_target,
110 const int compIdx)
const;
115 const int compIdx)
const;
119 const int seg_upwind,
120 const size_t comp_idx)
const;
130 const int comp_idx)
const;
133 EvalWell
getQs(
const int comp_idx)
const;
139 const std::array<EvalWell,numWellEq>&
eval(
const int idx)
const
140 {
return evaluation_[idx]; }
144 void processFractions(
const int seg);
147 EvalWell volumeFraction(
const int seg,
148 const unsigned compIdx)
const;
152 std::vector<std::array<double, numWellEq>> value_;
156 std::vector<std::array<EvalWell, numWellEq>> evaluation_;
Definition: DeferredLogger.hpp:57
Definition: MultisegmentWellGeneric.hpp:42
Definition: MultisegmentWellPrimaryVariables.hpp:42
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
Definition: MultisegmentWellPrimaryVariables.cpp:603
void resize(const int numSegments)
Resize values and evaluations.
Definition: MultisegmentWellPrimaryVariables.cpp:46
void copyToWellState(const MultisegmentWellGeneric< Scalar > &mswell, const double rho, const bool stop_or_zero_rate_target, WellState &well_state, DeferredLogger &deferred_logger) const
Copy values to well state.
Definition: MultisegmentWellPrimaryVariables.cpp:208
void updateNewton(const BVectorWell &dwells, const double relaxation_factor, const double DFLimit, const bool stop_or_zero_rate_target, const double max_pressure_change)
Update values from newton update vector.
Definition: MultisegmentWellPrimaryVariables.cpp:155
EvalWell getSegmentPressure(const int seg) const
Get pressure for a segment.
Definition: MultisegmentWellPrimaryVariables.cpp:578
EvalWell getBhp() const
Get bottomhole pressure.
Definition: MultisegmentWellPrimaryVariables.cpp:586
EvalWell getWQTotal() const
Get WQTotal.
Definition: MultisegmentWellPrimaryVariables.cpp:611
EvalWell volumeFractionScaled(const int seg, const int compIdx) const
Returns scaled volume fraction for a component in a segment.
Definition: MultisegmentWellPrimaryVariables.cpp:505
void update(const WellState &well_state, const bool stop_or_zero_rate_target)
Copy values from well state.
Definition: MultisegmentWellPrimaryVariables.cpp:67
EvalWell getSegmentRateUpwinding(const int seg, const int seg_upwind, const size_t comp_idx) const
Returns upwinding rate for a component in a segment.
Definition: MultisegmentWellPrimaryVariables.cpp:538
void init()
Initialize evaluations from values.
Definition: MultisegmentWellPrimaryVariables.cpp:54
EvalWell surfaceVolumeFraction(const int seg, const int compIdx) const
Returns surface volume fraction for a component in a segment.
Definition: MultisegmentWellPrimaryVariables.cpp:522
const std::array< EvalWell, numWellEq > & eval(const int idx) const
Returns a const ref to an evaluation.
Definition: MultisegmentWellPrimaryVariables.hpp:139
EvalWell getSegmentRate(const int seg, const int comp_idx) const
Get rate for a component in a segment.
Definition: MultisegmentWellPrimaryVariables.cpp:594
Definition: WellInterfaceIndices.hpp:33
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition: WellState.hpp:60
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27