#include <WellContributions.hpp>

Inheritance diagram for Opm::WellContributions:
Inheritance graph

Public Types

enum class  MatrixType { C , D , B }
 StandardWell has C, D and B matrices that need to be copied. More...
 
using UMFPackIndex = SuiteSparse_long
 

Public Member Functions

unsigned int getNumWells ()
 
void addNumBlocks (unsigned int numBlocks)
 
void alloc ()
 Allocate memory for the StandardWells. More...
 
virtual ~WellContributions ()
 Empty destructor. More...
 
void setBlockSize (unsigned int dim, unsigned int dim_wells)
 
void setVectorSize (unsigned N)
 
void addMatrix (MatrixType type, int *colIndices, double *values, unsigned int val_size)
 
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)
 

Static Public Member Functions

static std::unique_ptr< WellContributionscreate (const std::string &accelerator_mode, bool useWellConn)
 

Protected Member Functions

virtual void APIalloc ()
 API specific allocation. More...
 
virtual void APIaddMatrix (MatrixType, int *, double *, unsigned int)
 Api specific upload of matrix. More...
 

Protected Attributes

bool allocated = false
 
unsigned int N
 
unsigned int dim
 
unsigned int dim_wells
 
unsigned int num_blocks = 0
 
unsigned int num_std_wells = 0
 
unsigned int num_ms_wells = 0
 
unsigned int num_blocks_so_far = 0
 
unsigned int num_std_wells_so_far = 0
 
std::vector< unsigned int > val_pointers
 
std::vector< std::unique_ptr< MultisegmentWellContribution > > multisegments
 

Detailed Description

This class serves to eliminate the need to include the WellContributions into the matrix (with –matrix-add-well-contributions=true) for the cusparseSolver or openclSolver. If the –matrix-add-well-contributions commandline parameter is true, this class should still be used, but be empty. StandardWell and MultisegmentWell are supported for both cusparseSolver and openclSolver. A single instance (or pointer) of this class is passed to the BdaSolver. For StandardWell, this class contains all the data and handles the computation. For MultisegmentWell, the vector 'multisegments' contains all the data. For more information, check the MultisegmentWellContribution class. A StandardWell uses C, D and B and performs y -= (C^T * (D^-1 * (B*x))) B and C are vectors, disguised as matrices and contain blocks of StandardWell::numEq by StandardWell::numStaticWellEq D is a block, disguised as matrix, the square block has size StandardWell::numStaticWellEq. D is actually stored as D^-1 B*x and D*B*x are a vector with numStaticWellEq doubles C*D*B*x is a blocked matrix with a symmetric sparsity pattern, contains square blocks with size numEq. For every columnindex i, j in StandardWell::duneB_, there is a block on (i, j) in C*D*B*x.

This class is used in 3 phases:

  • get total size of all wellcontributions that must be stored here
  • allocate memory
  • copy data of wellcontributions

Member Typedef Documentation

◆ UMFPackIndex

using Opm::WellContributions::UMFPackIndex = SuiteSparse_long

Member Enumeration Documentation

◆ MatrixType

StandardWell has C, D and B matrices that need to be copied.

Enumerator

Constructor & Destructor Documentation

◆ ~WellContributions()

virtual Opm::WellContributions::~WellContributions ( )
virtual

Empty destructor.

Member Function Documentation

◆ addMatrix()

void Opm::WellContributions::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

Parameters
[in]typeindicate if C, D or B is sent
[in]colIndicescolumnindices of blocks in C or B, ignored for D
[in]valuesarray of nonzeroes
[in]val_sizenumber of blocks in C or B, ignored for D

◆ addMultisegmentWellContribution()

void Opm::WellContributions::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 destructor Matrices C and B are passed in Blocked CSR, matrix D in CSC

Parameters
[in]dimsize of blocks in vectors x and y, equal to MultisegmentWell::numEq
[in]dim_wellssize of blocks of C, B and D, equal to MultisegmentWell::numWellEq
[in]Mbnumber of blockrows in C, B and D
[in]Bvaluesnonzero values of matrix B
[in]BcolIndicescolumnindices of blocks of matrix B
[in]BrowPointersrowpointers of matrix B
[in]DnumBlocksnumber of blocks in D
[in]Dvaluesnonzero values of matrix D
[in]DcolPointerscolumnpointers of matrix D
[in]DrowIndicesrowindices of matrix D
[in]Cvaluesnonzero values of matrix C

◆ addNumBlocks()

void Opm::WellContributions::addNumBlocks ( unsigned int  numBlocks)

Indicate how large the next StandardWell is, this function cannot be called after alloc() is called

Parameters
[in]numBlocksnumber of blocks in C and B of next StandardWell

Referenced by Opm::BlackoilWellModel< TypeTag >::getWellContributions().

◆ alloc()

void Opm::WellContributions::alloc ( )

Allocate memory for the StandardWells.

Referenced by Opm::BlackoilWellModel< TypeTag >::getWellContributions().

◆ APIaddMatrix()

virtual void Opm::WellContributions::APIaddMatrix ( MatrixType  ,
int *  ,
double *  ,
unsigned int   
)
inlineprotectedvirtual

Api specific upload of matrix.

Reimplemented in Opm::WellContributionsCuda, Opm::WellContributionsOCL, and Opm::WellContributionsRocsparse.

◆ APIalloc()

virtual void Opm::WellContributions::APIalloc ( )
inlineprotectedvirtual

◆ create()

static std::unique_ptr< WellContributions > Opm::WellContributions::create ( const std::string &  accelerator_mode,
bool  useWellConn 
)
static

◆ getNumWells()

unsigned int Opm::WellContributions::getNumWells ( )
inline

References num_ms_wells, and num_std_wells.

◆ setBlockSize()

void Opm::WellContributions::setBlockSize ( unsigned int  dim,
unsigned int  dim_wells 
)

Indicate how large the blocks of the StandardWell (C and B) are

Parameters
[in]dimnumber of columns
[in]dim_wellsnumber of rows

Referenced by Opm::BlackoilWellModel< TypeTag >::getWellContributions().

◆ setVectorSize()

void Opm::WellContributions::setVectorSize ( unsigned  N)

Set size of vector that the wells are applied to

Parameters
[in]Nsize of vector

Member Data Documentation

◆ allocated

bool Opm::WellContributions::allocated = false
protected

◆ dim

unsigned int Opm::WellContributions::dim
protected

◆ dim_wells

unsigned int Opm::WellContributions::dim_wells
protected

◆ multisegments

std::vector<std::unique_ptr<MultisegmentWellContribution> > Opm::WellContributions::multisegments
protected

◆ N

unsigned int Opm::WellContributions::N
protected

◆ num_blocks

unsigned int Opm::WellContributions::num_blocks = 0
protected

◆ num_blocks_so_far

unsigned int Opm::WellContributions::num_blocks_so_far = 0
protected

◆ num_ms_wells

unsigned int Opm::WellContributions::num_ms_wells = 0
protected

Referenced by getNumWells().

◆ num_std_wells

unsigned int Opm::WellContributions::num_std_wells = 0
protected

Referenced by getNumWells().

◆ num_std_wells_so_far

unsigned int Opm::WellContributions::num_std_wells_so_far = 0
protected

◆ val_pointers

std::vector<unsigned int> Opm::WellContributions::val_pointers
protected

The documentation for this class was generated from the following file: