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;
42class WellContributions;
44class WellInterfaceGeneric;
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>;
80 const std::vector<int>& cells,
88 void apply(
const BVector& x, BVector&
Ax)
const;
91 void apply(BVector&
r)
const;
97 BVectorWell
solve()
const;
103#if COMPILE_BDA_BRIDGE
109 template<
class SparseMatrixAdapter>
113 template<
class PressureMatrix>
115 const BVector& weights,
116 const int pressureVarIndex,
131 OffDiagMatWell duneB_;
132 OffDiagMatWell duneC_;
139 mutable std::shared_ptr<Dune::UMFPack<DiagMatWell>> duneDSolver_;
142 BVectorWell resWell_;
Definition AquiferInterface.hpp:35
Class administering assembler access to equation system.
Definition MultisegmentWellAssemble.cpp:45
Definition MultisegmentWellEquations.hpp:49
void clear()
Set all coefficients to 0.
Definition MultisegmentWellEquations.cpp:128
void extractCPRPressureMatrix(PressureMatrix &jacobian, const BVector &weights, const int pressureVarIndex, const bool, const WellInterfaceGeneric &well, const int seg_pressure_var_ind, const WellState &well_state) const
Extract CPR pressure matrix.
Definition MultisegmentWellEquations.cpp:303
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
Definition MultisegmentWellEquations.cpp:139
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
Definition MultisegmentWellEquations.cpp:267
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Definition MultisegmentWellEquations.cpp:186
void createSolver()
Compute the LU-decomposition of D matrix.
Definition MultisegmentWellEquations.cpp:163
BVectorWell solve() const
Apply inverted D matrix to residual and return result.
Definition MultisegmentWellEquations.cpp:179
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition MultisegmentWellEquations.hpp:123
void init(const int num_cells, 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:56
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:52
Definition WellInterfaceGeneric.hpp:51
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