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>
29#include <opm/input/eclipse/Units/Units.hpp>
39template<
class Scalar>
class MultisegmentWellGeneric;
40template<
class Flu
idSystem,
class Indices>
class WellInterfaceIndices;
41template<
class Scalar>
class WellState;
43template<
class Flu
idSystem,
class Indices>
60 static constexpr bool has_wfrac_variable = Indices::waterEnabled && Indices::oilEnabled;
61 static constexpr bool has_gfrac_variable = Indices::gasEnabled && Indices::numPhases > 1;
63 static constexpr int WQTotal = 0;
64 static constexpr int WFrac = has_wfrac_variable ? 1 : -1000;
65 static constexpr int GFrac = has_gfrac_variable ? has_wfrac_variable + 1 : -1000;
66 static constexpr int SPres = has_wfrac_variable + has_gfrac_variable + 1;
69 static constexpr int numWellEq = Indices::numPhases + 1;
71 using Scalar =
typename FluidSystem::Scalar;
72 using EvalWell = DenseAd::Evaluation<Scalar, Indices::numEq + numWellEq>;
75 using BVectorWell =
typename Equations::BVectorWell;
134 const std::array<EvalWell, numWellEq>&
eval(
const int idx)
const
135 {
return evaluation_[
idx]; }
138 const std::array<Scalar, numWellEq>&
value(
const int idx)
const
139 {
return value_[
idx]; }
150 void setEvaluationsFromValues();
153 void processFractions(
const int seg);
156 EvalWell volumeFraction(
const int seg,
161 std::vector<std::array<Scalar, numWellEq>> value_;
165 std::vector<std::array<EvalWell, numWellEq>> evaluation_;
169 static constexpr double bhp_lower_limit = 1. * unit::barsa - 1. * unit::Pascal;
170 static constexpr double seg_pres_lower_limit = 0.;
Definition DeferredLogger.hpp:57
Definition MultisegmentWellGeneric.hpp:39
Definition MultisegmentWellPrimaryVariables.hpp:45
void resize(const int numSegments)
Resize values and evaluations.
Definition MultisegmentWellPrimaryVariables.cpp:49
EvalWell volumeFractionScaled(const int seg, const int compIdx) const
Returns scaled volume fraction for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:554
void outputLowLimitPressureSegments(DeferredLogger &deferred_logger) const
output the segments with pressure close to lower pressure limit for debugging purpose
Definition MultisegmentWellPrimaryVariables.cpp:668
EvalWell getSegmentRate(const int seg, const int comp_idx) const
Get rate for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:643
EvalWell getSegmentRateUpwinding(const int seg, const int seg_upwind, const int comp_idx) const
Returns upwinding rate for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:587
EvalWell getWQTotal() const
Get WQTotal.
Definition MultisegmentWellPrimaryVariables.cpp:660
void copyToWellState(const MultisegmentWellGeneric< Scalar > &mswell, const Scalar rho, WellState< Scalar > &well_state, const SummaryState &summary_state, DeferredLogger &deferred_logger) const
Copy values to well state.
Definition MultisegmentWellPrimaryVariables.cpp:217
const std::array< EvalWell, numWellEq > & eval(const int idx) const
Returns a const ref to an array of evaluations.
Definition MultisegmentWellPrimaryVariables.hpp:134
EvalWell getBhp() const
Get bottomhole pressure.
Definition MultisegmentWellPrimaryVariables.cpp:635
EvalWell surfaceVolumeFraction(const int seg, const int compIdx) const
Returns surface volume fraction for a component in a segment.
Definition MultisegmentWellPrimaryVariables.cpp:571
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:142
EvalWell getQs(const int comp_idx) const
Returns scaled rate for a component.
Definition MultisegmentWellPrimaryVariables.cpp:652
void update(const WellState< Scalar > &well_state, const bool stop_or_zero_rate_target)
Copy values from well state.
Definition MultisegmentWellPrimaryVariables.cpp:70
EvalWell getSegmentPressure(const int seg) const
Get pressure for a segment.
Definition MultisegmentWellPrimaryVariables.cpp:627
void updateNewton(const BVectorWell &dwells, const Scalar relaxation_factor, const Scalar DFLimit, const bool stop_or_zero_rate_target, const Scalar max_pressure_change)
Update values from newton update vector.
Definition MultisegmentWellPrimaryVariables.cpp:160
const std::array< Scalar, numWellEq > & value(const int idx) const
Returns a value array.
Definition MultisegmentWellPrimaryVariables.hpp:138
Definition WellInterfaceIndices.hpp:34
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:66
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242