23#ifndef OPM_STANDARDWELL_EQUATIONS_HEADER_INCLUDED
24#define OPM_STANDARDWELL_EQUATIONS_HEADER_INCLUDED
26#include <opm/simulators/utils/ParallelCommunication.hpp>
27#include <opm/simulators/wells/WellHelpers.hpp>
28#include <opm/common/TimingMacros.hpp>
29#include <dune/common/dynmatrix.hh>
30#include <dune/common/dynvector.hh>
31#include <dune/istl/bcrsmatrix.hh>
32#include <dune/istl/bvector.hh>
37class ParallelWellInfo;
38template<
class Scalar,
int numEq>
class StandardWellEquationAccess;
40class WellContributions;
42class WellInterfaceGeneric;
45template<
class Scalar,
int numEq>
54 using VectorBlockWellType = Dune::DynamicVector<Scalar>;
55 using BVectorWell = Dune::BlockVector<VectorBlockWellType>;
58 using DiagMatrixBlockWellType = Dune::DynamicMatrix<Scalar>;
59 using DiagMatWell = Dune::BCRSMatrix<DiagMatrixBlockWellType>;
62 using OffDiagMatrixBlockWellType = Dune::DynamicMatrix<Scalar>;
63 using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
66 using BVector = Dune::BlockVector<Dune::FieldVector<Scalar,numEq>>;
78 const std::vector<int>& cells);
84 void apply(
const BVector& x, BVector&
Ax)
const;
87 void apply(BVector&
r)
const;
101 void extract(
const int numStaticWellEq,
106 template<
class SparseMatrixAdapter>
110 template<
class PressureMatrix>
112 const BVector& weights,
113 const int pressureVarIndex,
135 OffDiagMatWell duneB_;
136 OffDiagMatWell duneC_;
138 DiagMatWell invDuneD_;
145 BVectorWell resWell_;
148 mutable BVectorWell Bx_;
149 mutable BVectorWell invDrw_;
Definition AquiferInterface.hpp:35
Class encapsulating some information about parallel wells.
Definition ParallelWellInfo.hpp:184
Class administering assembler access to equation system.
Definition StandardWellAssemble.cpp:45
Definition StandardWellEquations.hpp:47
const BVectorWell & residual() const
Returns a const reference to the residual.
Definition StandardWellEquations.hpp:126
void extract(SparseMatrixAdapter &jacobian) const
Add the matrices of this well to the sparse matrix adapter.
Definition StandardWellEquations.cpp:251
void invert()
Invert D matrix.
Definition StandardWellEquations.cpp:160
void apply(const BVector &x, BVector &Ax) const
Apply linear operator to vector.
Definition StandardWellEquations.cpp:130
void init(const int num_cells, const int numWellEq, const int numPerfs, const std::vector< int > &cells)
Setup sparsity pattern for the matrices.
Definition StandardWellEquations.cpp:55
unsigned int getNumBlocks() const
Get the number of blocks of the C and B matrices.
Definition StandardWellEquations.cpp:277
void sumDistributed(Parallel::Communication comm)
Sum with off-process contribution.
Definition StandardWellEquations.cpp:396
void solve(BVectorWell &dx_well) const
Apply inverted D matrix to residual and store in vector.
Definition StandardWellEquations.cpp:175
void recoverSolutionWell(const BVector &x, BVectorWell &xw) const
Recover well solution.
Definition StandardWellEquations.cpp:183
void clear()
Set all coefficients to 0.
Definition StandardWellEquations.cpp:121
void extractCPRPressureMatrix(PressureMatrix &jacobian, const BVector &weights, const int pressureVarIndex, const bool use_well_weights, const WellInterfaceGeneric &well, const int bhp_var_index, const WellState &well_state) const
Extract CPR pressure matrix.
Definition StandardWellEquations.cpp:285
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
A wrapper around the B matrix for distributed wells.
Definition WellHelpers.hpp:53
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27