#include <WellContributions.hpp>
|
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) |
|
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
◆ UMFPackIndex
◆ MatrixType
StandardWell has C, D and B matrices that need to be copied.
◆ ~WellContributions()
virtual Opm::WellContributions::~WellContributions |
( |
| ) |
|
|
virtual |
◆ 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] | type | indicate if C, D or B is sent |
[in] | colIndices | columnindices of blocks in C or B, ignored for D |
[in] | values | array of nonzeroes |
[in] | val_size | number 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] | dim | size of blocks in vectors x and y, equal to MultisegmentWell::numEq |
[in] | dim_wells | size of blocks of C, B and D, equal to MultisegmentWell::numWellEq |
[in] | Mb | number of blockrows in C, B and D |
[in] | Bvalues | nonzero values of matrix B |
[in] | BcolIndices | columnindices of blocks of matrix B |
[in] | BrowPointers | rowpointers of matrix B |
[in] | DnumBlocks | number of blocks in D |
[in] | Dvalues | nonzero values of matrix D |
[in] | DcolPointers | columnpointers of matrix D |
[in] | DrowIndices | rowindices of matrix D |
[in] | Cvalues | nonzero values of matrix C |
◆ addNumBlocks()
void Opm::WellContributions::addNumBlocks |
( |
unsigned int |
numBlocks | ) |
|
◆ alloc()
void Opm::WellContributions::alloc |
( |
| ) |
|
◆ APIaddMatrix()
virtual void Opm::WellContributions::APIaddMatrix |
( |
MatrixType |
, |
|
|
int * |
, |
|
|
double * |
, |
|
|
unsigned int |
|
|
) |
| |
|
inlineprotectedvirtual |
◆ 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 |
◆ setBlockSize()
void Opm::WellContributions::setBlockSize |
( |
unsigned int |
dim, |
|
|
unsigned int |
dim_wells |
|
) |
| |
◆ setVectorSize()
void Opm::WellContributions::setVectorSize |
( |
unsigned |
N | ) |
|
Set size of vector that the wells are applied to - Parameters
-
◆ allocated
bool Opm::WellContributions::allocated = false |
|
protected |
◆ dim
unsigned int Opm::WellContributions::dim |
|
protected |
◆ dim_wells
unsigned int Opm::WellContributions::dim_wells |
|
protected |
◆ multisegments
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 |
◆ num_std_wells
unsigned int Opm::WellContributions::num_std_wells = 0 |
|
protected |
◆ 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:
|