23#ifndef OPM_WELLHELPERS_HEADER_INCLUDED
24#define OPM_WELLHELPERS_HEADER_INCLUDED
26#include <dune/istl/bcrsmatrix.hh>
27#include <dune/common/dynmatrix.hh>
31template<
class Scalar>
class ParallelWellInfo;
32struct WellProductionControls;
33struct WellInjectionControls;
34enum class WellProducerCMode;
35enum class WellInjectorCMode;
37namespace wellhelpers {
49template<
typename Scalar>
53 using Block = Dune::DynamicMatrix<Scalar>;
54 using Matrix = Dune::BCRSMatrix<Block>;
60 template<
class X,
class Y>
61 void mv (
const X& x, Y& y)
const;
64 template<
class X,
class Y>
65 void mmv (
const X& x, Y& y)
const;
73Scalar computeHydrostaticCorrection(
const Scalar well_ref_depth,
74 const Scalar vfp_ref_depth,
76 const Scalar gravity);
79template<
typename Scalar,
typename Comm>
80void sumDistributedWellEntries(Dune::DynamicMatrix<Scalar>& mat,
81 Dune::DynamicVector<Scalar>& vec,
87template <
class DenseMatrix>
88DenseMatrix transposeDenseDynMatrix(
const DenseMatrix& M);
91bool rateControlWithZeroProdTarget(
const WellProductionControls& controls,
92 WellProducerCMode mode);
95bool rateControlWithZeroInjTarget(
const WellInjectionControls& controls,
96 WellInjectorCMode mode);
Class encapsulating some information about parallel wells.
Definition WGState.hpp:31
A wrapper around the B matrix for distributed wells.
Definition WellHelpers.hpp:51
void mv(const X &x, Y &y) const
y = A x
Definition WellHelpers.cpp:52
void mmv(const X &x, Y &y) const
y = A x
Definition WellHelpers.cpp:109
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37