20#ifndef WELLCONTRIBUTIONS_HEADER_INCLUDED
21#define WELLCONTRIBUTIONS_HEADER_INCLUDED
26#if HAVE_SUITESPARSE_UMFPACK
29#include <dune/common/version.hh>
33class MultisegmentWellContribution;
54 static std::unique_ptr<WellContributions> create(
const std::string& accelerator_mode,
bool useWellConn);
56 using UMFPackIndex = SuiteSparse_long;
65 bool allocated =
false;
69 unsigned int dim_wells;
70 unsigned int num_blocks = 0;
71 unsigned int num_std_wells = 0;
72 unsigned int num_ms_wells = 0;
73 unsigned int num_blocks_so_far = 0;
74 unsigned int num_std_wells_so_far = 0;
75 std::vector<unsigned int> val_pointers;
77 std::vector<std::unique_ptr<MultisegmentWellContribution>> multisegments;
80 unsigned int getNumWells(){
81 return num_std_wells + num_ms_wells;
97 void setBlockSize(
unsigned int dim,
unsigned int dim_wells);
125 std::vector<double> &Bvalues, std::vector<unsigned int> &BcolIndices, std::vector<unsigned int> &BrowPointers,
126 unsigned int DnumBlocks,
double *Dvalues,
127 UMFPackIndex *DcolPointers, UMFPackIndex *DrowIndices,
128 std::vector<double> &Cvalues);
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:52
virtual ~WellContributions()
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:110
virtual void APIalloc()
API specific allocation.
Definition WellContributions.hpp:131
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:89
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:147
MatrixType
StandardWell has C, D and B matrices that need to be copied.
Definition WellContributions.hpp:58
void alloc()
Allocate memory for the StandardWells.
Definition WellContributions.cpp:137
void setVectorSize(unsigned N)
Set size of vector that the wells are applied to.
Definition WellContributions.cpp:124
void addNumBlocks(unsigned int numBlocks)
Indicate how large the next StandardWell is, this function cannot be called after alloc() is called.
Definition WellContributions.cpp:128
virtual void APIaddMatrix(MatrixType, int *, double *, unsigned int)
Api specific upload of matrix.
Definition WellContributions.hpp:134
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27