WellContributions.hpp
Go to the documentation of this file.
1/*
2 Copyright 2020 Equinor ASA
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef WELLCONTRIBUTIONS_HEADER_INCLUDED
21#define WELLCONTRIBUTIONS_HEADER_INCLUDED
22
23#include <memory>
24#include <vector>
25
26#if HAVE_SUITESPARSE_UMFPACK
27#include<umfpack.h>
28#endif
29#include <dune/common/version.hh>
30
31namespace Opm {
32
33template<class Scalar> class MultisegmentWellContribution;
34
40
51template<class Scalar>
53{
54public:
55 static std::unique_ptr<WellContributions> create(const std::string& accelerator_mode, bool useWellConn);
56
57 using UMFPackIndex = SuiteSparse_long;
59 enum class MatrixType {
60 C,
61 D,
62 B
63 };
64
65protected:
66 bool allocated = false;
67
68 unsigned int N; // number of rows (not blockrows) in vectors x and y
69 unsigned int dim; // number of columns in blocks in B and C, equal to StandardWell::numEq
70 unsigned int dim_wells; // number of rows in blocks in B and C, equal to StandardWell::numStaticWellEq
71 unsigned int num_blocks = 0; // total number of blocks in all wells
72 unsigned int num_std_wells = 0; // number of StandardWells in this object
73 unsigned int num_ms_wells = 0; // number of MultisegmentWells in this object, must equal multisegments.size()
74 unsigned int num_blocks_so_far = 0; // keep track of where next data is written
75 unsigned int num_std_wells_so_far = 0; // keep track of where next data is written
76 std::vector<unsigned int> val_pointers; // val_pointers[wellID] == index of first block for this well in Ccols and Bcols
77
78 std::vector<std::unique_ptr<MultisegmentWellContribution<Scalar>>> multisegments;
79
80public:
81 unsigned int getNumWells(){
83 }
84
87 void addNumBlocks(unsigned int numBlocks);
88
90 void alloc();
91
94
98 void setBlockSize(unsigned int dim, unsigned int dim_wells);
99
102 void setVectorSize(unsigned N);
103
109 void addMatrix(MatrixType type, int* colIndices, Scalar* values, unsigned int val_size);
110
125 unsigned int dim_wells,
126 unsigned int Mb,
127 std::vector<Scalar>& Bvalues,
128 std::vector<unsigned int>& BcolIndices,
129 std::vector<unsigned int>& BrowPointers,
130 unsigned int DnumBlocks,
131 Scalar* Dvalues,
132 UMFPackIndex* DcolPointers,
133 UMFPackIndex* DrowIndices,
134 std::vector<Scalar>& Cvalues);
135protected:
137 virtual void APIalloc() {}
138
140 virtual void APIaddMatrix(MatrixType, int*, Scalar*, unsigned int) {}
141};
142
143} //namespace Opm
144
145#endif
Definition: WellContributions.hpp:53
void addMatrix(MatrixType type, int *colIndices, Scalar *values, unsigned int val_size)
void addMultisegmentWellContribution(unsigned int dim, unsigned int dim_wells, unsigned int Mb, std::vector< Scalar > &Bvalues, std::vector< unsigned int > &BcolIndices, std::vector< unsigned int > &BrowPointers, unsigned int DnumBlocks, Scalar *Dvalues, UMFPackIndex *DcolPointers, UMFPackIndex *DrowIndices, std::vector< Scalar > &Cvalues)
void alloc()
Allocate memory for the StandardWells.
void setBlockSize(unsigned int dim, unsigned int dim_wells)
unsigned int dim_wells
Definition: WellContributions.hpp:70
void addNumBlocks(unsigned int numBlocks)
SuiteSparse_long UMFPackIndex
Definition: WellContributions.hpp:57
unsigned int num_std_wells_so_far
Definition: WellContributions.hpp:75
void setVectorSize(unsigned N)
unsigned int N
Definition: WellContributions.hpp:68
virtual void APIaddMatrix(MatrixType, int *, Scalar *, unsigned int)
Api specific upload of matrix.
Definition: WellContributions.hpp:140
unsigned int num_ms_wells
Definition: WellContributions.hpp:73
unsigned int num_std_wells
Definition: WellContributions.hpp:72
unsigned int num_blocks
Definition: WellContributions.hpp:71
bool allocated
Definition: WellContributions.hpp:66
MatrixType
StandardWell has C, D and B matrices that need to be copied.
Definition: WellContributions.hpp:59
std::vector< std::unique_ptr< MultisegmentWellContribution< Scalar > > > multisegments
Definition: WellContributions.hpp:78
unsigned int num_blocks_so_far
Definition: WellContributions.hpp:74
static std::unique_ptr< WellContributions > create(const std::string &accelerator_mode, bool useWellConn)
unsigned int getNumWells()
Definition: WellContributions.hpp:81
virtual ~WellContributions()
Empty destructor.
unsigned int dim
Definition: WellContributions.hpp:69
virtual void APIalloc()
API specific allocation.
Definition: WellContributions.hpp:137
std::vector< unsigned int > val_pointers
Definition: WellContributions.hpp:76
Definition: blackoilboundaryratevector.hh:37