22#ifndef OPM_MULTISEGMENTWELL_EQUATIONS_HEADER_INCLUDED
23#define OPM_MULTISEGMENTWELL_EQUATIONS_HEADER_INCLUDED
25#include <opm/simulators/utils/ParallelCommunication.hpp>
26#include <opm/simulators/wells/ParallelWellInfo.hpp>
27#include <opm/simulators/wells/MSWellHelpers.hpp>
28#include <dune/common/fmatrix.hh>
29#include <dune/common/fvector.hh>
30#include <dune/istl/bcrsmatrix.hh>
31#include <dune/istl/bvector.hh>
48template<
typename Scalar,
typename IndexTraits>
class WellState;
50template<
class Scalar,
typename IndexTraits,
int numWellEq,
int numEq>
51class MultisegmentWellEquations
59 using VectorBlockWellType = Dune::FieldVector<Scalar,numWellEq>;
60 using BVectorWell = Dune::BlockVector<VectorBlockWellType>;
62 using VectorBlockType = Dune::FieldVector<Scalar,numEq>;
63 using BVector = Dune::BlockVector<VectorBlockType>;
66 using DiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numWellEq>;
67 using DiagMatWell = Dune::BCRSMatrix<DiagMatrixBlockWellType>;
70 using OffDiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numEq>;
71 using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
80 void init(
const int numPerfs,
81 const std::vector<int>& cells,
82 const std::vector<std::vector<int>>& segment_inlets,
83 const std::vector<std::vector<int>>& segment_perforations);
89 void apply(
const BVector& x, BVector& Ax)
const;
101 BVectorWell
solve(
const BVectorWell& rhs)
const;
107#if COMPILE_GPU_BRIDGE
113 template<
class SparseMatrixAdapter>
114 void extract(SparseMatrixAdapter& jacobian)
const;
117 template<
class PressureMatrix>
119 const BVector& weights,
120 const int pressureVarIndex,
123 const int seg_pressure_var_ind,
138 OffDiagMatWell duneB_;
139 OffDiagMatWell duneC_;
146 mutable std::shared_ptr<Dune::UMFPack<DiagMatWell>> duneDSolver_;
149 BVectorWell resWell_;
154 std::vector<int> cells_;
Definition MSWellHelpers.hpp:29
Class administering assembler access to equation system.
Definition MultisegmentWellAssemble.cpp:44
void sumDistributed(Parallel::Communication comm)
Sum with off-process contribution.
Definition MultisegmentWellEquations.cpp:441
BVectorWell solve(const BVectorWell &rhs) const
Apply inverted D matrix to rhs and return result.
Definition MultisegmentWellEquations.cpp:220
BVectorWell solve() const
Apply inverted D matrix to residual and return result.
Definition MultisegmentWellEquations.cpp:210
void apply(BVector &r) const
Apply linear operator to vector.
Definition MultisegmentWellEquations.cpp:173
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition MultisegmentWellEquations.hpp:130
void extractCPRPressureMatrix(PressureMatrix &jacobian, const BVector &weights, const int pressureVarIndex, const bool, const WellInterfaceGeneric< Scalar, IndexTraits > &well, const int seg_pressure_var_ind, const WellState< Scalar, IndexTraits > &well_state) const
Extract CPR pressure matrix.
Definition MultisegmentWellEquations.cpp:356
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
Definition MultisegmentWellEquations.cpp:319
void clear()
Set all coefficients to 0.
Definition MultisegmentWellEquations.cpp:141
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
Definition MultisegmentWellEquations.cpp:152
void createSolver()
Compute the LU-decomposition of D matrix.
Definition MultisegmentWellEquations.cpp:189
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:63
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Definition MultisegmentWellEquations.cpp:230
Definition MultisegmentWellGeneric.hpp:39
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:198
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51
Definition WellInterfaceGeneric.hpp:53
The state of a set of wells, tailored for use by the fully implicit blackoil simulator.
Definition WellState.hpp:66
A wrapper around the B matrix for distributed MS wells.
Definition MSWellHelpers.hpp:54
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43