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>
49template<
typename Scalar,
typename IndexTraits>
class WellState;
51template<
class Scalar,
typename IndexTraits,
int numWellEq,
int numEq>
52class MultisegmentWellEquations
60 using VectorBlockWellType = Dune::FieldVector<Scalar,numWellEq>;
61 using BVectorWell = Dune::BlockVector<VectorBlockWellType>;
63 using VectorBlockType = Dune::FieldVector<Scalar,numEq>;
64 using BVector = Dune::BlockVector<VectorBlockType>;
67 using DiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numWellEq>;
68 using DiagMatWell = Dune::BCRSMatrix<DiagMatrixBlockWellType>;
71 using OffDiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numEq>;
72 using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
82 void init(
const int numPerfs,
83 const std::vector<int>& cells,
84 const std::vector<std::vector<int>>& segment_inlets,
85 const std::vector<std::vector<int>>& segment_perforations,
92 void apply(
const BVector& x, BVector& Ax)
const;
104 BVectorWell
solve(
const BVectorWell& rhs)
const;
110#if COMPILE_GPU_BRIDGE
116 template<
class SparseMatrixAdapter>
117 void extract(SparseMatrixAdapter& jacobian)
const;
120 template<
class PressureMatrix>
122 const BVector& weights,
123 const int pressureVarIndex,
126 const int seg_pressure_var_ind,
141 OffDiagMatWell duneB_;
142 OffDiagMatWell duneC_;
150 using UMFPackSolver = std::conditional_t<std::is_same_v<Scalar,double>,
151 std::shared_ptr<Dune::UMFPack<DiagMatWell>>,
153 mutable UMFPackSolver duneDSolver_;
156 BVectorWell resWell_;
158 const MultisegmentWellGeneric<Scalar, IndexTraits>& well_;
161 std::vector<int> cells_;
164 mswellhelpers::ParallellMSWellB<OffDiagMatWell> parallelB_;
Definition MSWellHelpers.hpp:29
Class administering assembler access to equation system.
Definition MultisegmentWellAssemble.cpp:45
void sumDistributed(Parallel::Communication comm)
Sum with off-process contribution.
Definition MultisegmentWellEquations.cpp:457
BVectorWell solve(const BVectorWell &rhs) const
Apply inverted D matrix to rhs and return result.
Definition MultisegmentWellEquations.cpp:230
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, const ParallelWellInfo< Scalar > ¶llel_well_info)
Setup sparsity pattern for the matrices.
Definition MultisegmentWellEquations.cpp:62
BVectorWell solve() const
Apply inverted D matrix to residual and return result.
Definition MultisegmentWellEquations.cpp:215
void apply(BVector &r) const
Apply linear operator to vector.
Definition MultisegmentWellEquations.cpp:177
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition MultisegmentWellEquations.hpp:133
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:372
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
Definition MultisegmentWellEquations.cpp:333
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:154
void createSolver()
Compute the LU-decomposition of D matrix.
Definition MultisegmentWellEquations.cpp:195
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Definition MultisegmentWellEquations.cpp:245
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
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45