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>
28#include <opm/input/eclipse/Schedule/SummaryState.hpp>
38template<
class Scalar>
class MultisegmentWellGeneric;
39template<
class Flu
idSystem,
class Indices,
class Scalar>
class WellInterfaceIndices;
42template<
class Flu
idSystem,
class Indices,
class Scalar>
59 static constexpr bool has_water = (Indices::waterSwitchIdx >= 0);
60 static constexpr bool has_gas = (Indices::compositionSwitchIdx >= 0);
61 static constexpr bool has_oil = (Indices::numPhases - has_gas - has_water) > 0;
66 static constexpr bool has_wfrac_variable = has_water && Indices::numPhases > 1;
67 static constexpr bool has_gfrac_variable = has_gas && has_oil;
69 static constexpr int WQTotal = 0;
70 static constexpr int WFrac = has_wfrac_variable ? 1 : -1000;
71 static constexpr int GFrac = has_gfrac_variable ? has_wfrac_variable + 1 : -1000;
72 static constexpr int SPres = has_wfrac_variable + has_gfrac_variable + 1;
75 static constexpr int numWellEq = Indices::numPhases + 1;
77 using EvalWell = DenseAd::Evaluation<
double, Indices::numEq + numWellEq>;
80 using BVectorWell =
typename Equations::BVectorWell;
142 const std::array<EvalWell, numWellEq>&
eval(
const int idx)
const
143 {
return evaluation_[
idx]; }
146 const std::array<Scalar, numWellEq>&
value(
const int idx)
const
147 {
return value_[
idx]; }
155 void processFractions(
const int seg);
158 EvalWell volumeFraction(
const int seg,
163 std::vector<std::array<double, numWellEq>> value_;
167 std::vector<std::array<EvalWell, numWellEq>> evaluation_;
Definition AquiferInterface.hpp:35
Definition DeferredLogger.hpp:57
Definition MultisegmentWellPrimaryVariables.hpp:44
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
Definition MultisegmentWellPrimaryVariables.cpp:618
void resize(const int numSegments)
Resize values and evaluations.
Definition MultisegmentWellPrimaryVariables.cpp:47
void setValue(const int idx, const std::array< Scalar, numWellEq > &val)
Set a value array. Note that this does not also set the corresponding evaluation.
Definition MultisegmentWellPrimaryVariables.hpp:150
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:156
EvalWell getSegmentPressure(const int seg) const
Get pressure for a segment.
Definition MultisegmentWellPrimaryVariables.cpp:593
EvalWell getBhp() const
Get bottomhole pressure.
Definition MultisegmentWellPrimaryVariables.cpp:601
EvalWell getWQTotal() const
Get WQTotal.
Definition MultisegmentWellPrimaryVariables.cpp:626
EvalWell volumeFractionScaled(const int seg, const int compIdx) const
Returns scaled volume fraction for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:520
const std::array< Scalar, numWellEq > & value(const int idx) const
Returns a value array.
Definition MultisegmentWellPrimaryVariables.hpp:146
void update(const WellState &well_state, const bool stop_or_zero_rate_target)
Copy values from well state.
Definition MultisegmentWellPrimaryVariables.cpp:68
void copyToWellState(const MultisegmentWellGeneric< Scalar > &mswell, const double rho, const bool stop_or_zero_rate_target, WellState &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Copy values to well state.
Definition MultisegmentWellPrimaryVariables.cpp:214
void init()
Initialize evaluations from values.
Definition MultisegmentWellPrimaryVariables.cpp:55
EvalWell surfaceVolumeFraction(const int seg, const int compIdx) const
Returns surface volume fraction for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:537
const std::array< EvalWell, numWellEq > & eval(const int idx) const
Returns a const ref to an array of evaluations.
Definition MultisegmentWellPrimaryVariables.hpp:142
EvalWell getSegmentRateUpwinding(const int seg, const int seg_upwind, const std::size_t comp_idx) const
Returns upwinding rate for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:553
EvalWell getSegmentRate(const int seg, const int comp_idx) const
Get rate for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:609
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