22#ifndef OPM_MULTISEGMENTWELL_EVAL_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_EVAL_HEADER_INCLUDED
25#include <opm/simulators/wells/MultisegmentWellEquations.hpp>
26#include <opm/simulators/wells/MultisegmentWellGeneric.hpp>
27#include <opm/simulators/wells/MultisegmentWellPrimaryVariables.hpp>
28#include <opm/simulators/wells/MultisegmentWellSegments.hpp>
29#include <opm/simulators/wells/ParallelWellInfo.hpp>
31#include <opm/material/densead/Evaluation.hpp>
38class ConvergenceReport;
42template<
class Flu
idSystem,
class Indices>
class WellInterfaceIndices;
43template<
typename Scalar,
typename IndexTraits>
class WellState;
45template<
typename Flu
idSystem,
typename Indices>
47 typename FluidSystem::IndexTraitsType>
50 using Scalar =
typename FluidSystem::Scalar;
51 using IndexTraits =
typename FluidSystem::IndexTraitsType;
53 static constexpr int numWellEq = PrimaryVariables::numWellEq;
54 static constexpr int SPres = PrimaryVariables::SPres;
55 static constexpr int WQTotal = PrimaryVariables::WQTotal;
60 using BVector =
typename Equations::BVector;
61 using BVectorWell =
typename Equations::BVectorWell;
66 using EvalWell =
typename PrimaryVariables::EvalWell;
67 using Eval = DenseAd::Evaluation<Scalar, Indices::numEq>;
78 void initMatrixAndVectors();
80 void assembleDefaultPressureEq(
const int seg,
85 void assembleICDPressureEq(
const int seg,
92 void assembleAccelerationAndHydroPressureLosses(
const int seg,
97 void assemblePressureEq(
const int seg,
106 const std::vector<Scalar>&
B_avg,
108 const Scalar max_residual_allowed,
116 std::pair<bool, std::vector<Scalar> >
117 getFiniteWellResiduals(
const std::vector<Scalar>&
B_avg,
131 void assembleAccelerationPressureLoss(
const int seg,
134 EvalWell pressureDropAutoICD(
const int seg,
138 EvalWell extendEval(
const Eval&
in)
const;
147 std::vector<Scalar> cell_perforation_depth_diffs_;
150 std::vector<Scalar> cell_perforation_pressure_diffs_;
Represents the convergence status of the whole simulator, to make it possible to query and store the ...
Definition ConvergenceReport.hpp:38
Definition DeferredLogger.hpp:57
Definition MultisegmentWellEval.hpp:48
Equations linSys_
The equation system.
Definition MultisegmentWellEval.hpp:142
const Equations & linSys() const
Returns a const reference to equation system.
Definition MultisegmentWellEval.hpp:71
PrimaryVariables primary_variables_
The primary variables.
Definition MultisegmentWellEval.hpp:143
MSWSegments segments_
Segment properties.
Definition MultisegmentWellEval.hpp:144
ConvergenceReport getWellConvergence(const WellState< Scalar, IndexTraits > &well_state, const std::vector< Scalar > &B_avg, DeferredLogger &deferred_logger, const Scalar max_residual_allowed, const Scalar tolerance_wells, const Scalar relaxed_inner_tolerance_flow_ms_well, const Scalar tolerance_pressure_ms_wells, const Scalar relaxed_inner_tolerance_pressure_ms_well, const bool relax_tolerance, const bool well_is_stopped) const
check whether the well equations get converged for this well
Definition MultisegmentWellEval.cpp:82
Definition MultisegmentWellGeneric.hpp:39
Definition MultisegmentWellPrimaryVariables.hpp:45
Definition MultisegmentWellSegments.hpp:45
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:198
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 blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240