22#ifndef OPM_MULTISEGMENTWELL_EQUATIONS_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_EQUATIONS_HEADER_INCLUDED
25#include <dune/common/fmatrix.hh>
26#include <dune/common/fvector.hh>
27#include <dune/istl/bcrsmatrix.hh>
28#include <dune/istl/bvector.hh>
33template<
class M>
class UMFPack;
39template<
class Scalar,
int numWellEq,
int numEq>
class MultisegmentWellEquationAccess;
40template<
class Scalar>
class MultisegmentWellGeneric;
42template<
class Scalar>
class WellContributions;
44template<
class Scalar>
class WellInterfaceGeneric;
45template<
class Scalar>
class WellState;
47template<
class Scalar,
int numWellEq,
int numEq>
56 using VectorBlockWellType = Dune::FieldVector<Scalar,numWellEq>;
57 using BVectorWell = Dune::BlockVector<VectorBlockWellType>;
59 using VectorBlockType = Dune::FieldVector<Scalar,numEq>;
60 using BVector = Dune::BlockVector<VectorBlockType>;
63 using DiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numWellEq>;
64 using DiagMatWell = Dune::BCRSMatrix<DiagMatrixBlockWellType>;
67 using OffDiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numEq>;
68 using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
77 void init(
const int numPerfs,
78 const std::vector<int>& cells,
86 void apply(
const BVector& x, BVector&
Ax)
const;
89 void apply(BVector&
r)
const;
95 BVectorWell
solve()
const;
98 BVectorWell
solve(
const BVectorWell& rhs)
const;
104#if COMPILE_BDA_BRIDGE
110 template<
class SparseMatrixAdapter>
111 void extract(SparseMatrixAdapter& jacobian)
const;
114 template<
class PressureMatrix>
116 const BVector& weights,
117 const int pressureVarIndex,
132 OffDiagMatWell duneB_;
133 OffDiagMatWell duneC_;
140 mutable std::shared_ptr<Dune::UMFPack<DiagMatWell>> duneDSolver_;
143 BVectorWell resWell_;
148 std::vector<int> cells_;
Class administering assembler access to equation system.
Definition MultisegmentWellAssemble.cpp:46
Definition MultisegmentWellEquations.hpp:49
void extractCPRPressureMatrix(PressureMatrix &jacobian, const BVector &weights, const int pressureVarIndex, const bool, const WellInterfaceGeneric< Scalar > &well, const int seg_pressure_var_ind, const WellState< Scalar > &well_state) const
Extract CPR pressure matrix.
Definition MultisegmentWellEquations.cpp:322
void clear()
Set all coefficients to 0.
Definition MultisegmentWellEquations.cpp:130
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
Definition MultisegmentWellEquations.cpp:141
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
Definition MultisegmentWellEquations.cpp:285
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Definition MultisegmentWellEquations.cpp:200
void createSolver()
Compute the LU-decomposition of D matrix.
Definition MultisegmentWellEquations.cpp:165
void init(const int numPerfs, const std::vector< int > &cells, const std::vector< std::vector< int > > &segment_inlets, const std::vector< std::vector< int > > &segment_perforations)
Setup sparsity pattern for the matrices.
Definition MultisegmentWellEquations.cpp:58
BVectorWell solve() const
Apply inverted D matrix to residual and return result.
Definition MultisegmentWellEquations.cpp:186
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition MultisegmentWellEquations.hpp:124
Definition MultisegmentWellGeneric.hpp:42
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:53
Definition WellInterfaceGeneric.hpp:52
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:62
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