20#ifndef WELLCONTRIBUTIONS_HEADER_INCLUDED
21#define WELLCONTRIBUTIONS_HEADER_INCLUDED
26#include <opm/simulators/linalg/bda/MultisegmentWellContribution.hpp>
27#if HAVE_SUITESPARSE_UMFPACK
30#include <dune/common/version.hh>
55 static std::unique_ptr<WellContributions> create(
const std::string& accelerator_mode,
bool useWellConn);
57#if DUNE_VERSION_NEWER(DUNE_ISTL, 2, 7)
58 using UMFPackIndex = SuiteSparse_long;
60 using UMFPackIndex = int;
70 bool allocated =
false;
74 unsigned int dim_wells;
75 unsigned int num_blocks = 0;
76 unsigned int num_std_wells = 0;
77 unsigned int num_ms_wells = 0;
78 unsigned int num_blocks_so_far = 0;
79 unsigned int num_std_wells_so_far = 0;
80 std::vector<unsigned int> val_pointers;
82 std::vector<std::unique_ptr<MultisegmentWellContribution>> multisegments;
85 unsigned int getNumWells(){
86 return num_std_wells + num_ms_wells;
102 void setBlockSize(
unsigned int dim,
unsigned int dim_wells);
130 std::vector<double> &Bvalues, std::vector<unsigned int> &BcolIndices, std::vector<unsigned int> &BrowPointers,
131 unsigned int DnumBlocks,
double *Dvalues,
132 UMFPackIndex *DcolPointers, UMFPackIndex *DrowIndices,
133 std::vector<double> &Cvalues);
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition: WellContributions.hpp:53
virtual ~WellContributions()=default
Empty destructor.
void setBlockSize(unsigned int dim, unsigned int dim_wells)
Indicate how large the blocks of the StandardWell (C and B) are.
Definition: WellContributions.cpp:92
virtual void APIalloc()
API specific allocation.
Definition: WellContributions.hpp:136
void addMatrix(MatrixType type, int *colIndices, double *values, unsigned int val_size)
Store a matrix in this object, in blocked csr format, can only be called after alloc() is called.
Definition: WellContributions.cpp:71
void addMultisegmentWellContribution(unsigned int dim, unsigned int dim_wells, unsigned int Mb, std::vector< double > &Bvalues, std::vector< unsigned int > &BcolIndices, std::vector< unsigned int > &BrowPointers, unsigned int DnumBlocks, double *Dvalues, UMFPackIndex *DcolPointers, UMFPackIndex *DrowIndices, std::vector< double > &Cvalues)
Add a MultisegmentWellContribution, actually creates an object on heap that is destroyed in the destr...
Definition: WellContributions.cpp:127
MatrixType
StandardWell has C, D and B matrices that need to be copied.
Definition: WellContributions.hpp:63
void alloc()
Allocate memory for the StandardWells.
Definition: WellContributions.cpp:117
void setVectorSize(unsigned N)
Set size of vector that the wells are applied to.
Definition: WellContributions.cpp:104
void addNumBlocks(unsigned int numBlocks)
Indicate how large the next StandardWell is, this function cannot be called after alloc() is called.
Definition: WellContributions.cpp:108
virtual void APIaddMatrix(MatrixType, int *, double *, unsigned int)
Api specific upload of matrix.
Definition: WellContributions.hpp:139
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27