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#include <umfpack.h>
27#include <dune/common/version.hh>
28
29namespace Opm {
30
31template<class Scalar> class MultisegmentWellContribution;
32
38
49template<class Scalar>
51{
52public:
53 static std::unique_ptr<WellContributions> create(const std::string& accelerator_mode, bool useWellConn);
54
55 using UMFPackIndex = SuiteSparse_long;
57 enum class MatrixType {
58 C,
59 D,
60 B
61 };
62
63protected:
64 bool allocated = false;
65
66 unsigned int N; // number of rows (not blockrows) in vectors x and y
67 unsigned int dim; // number of columns in blocks in B and C, equal to StandardWell::numEq
68 unsigned int dim_wells; // number of rows in blocks in B and C, equal to StandardWell::numStaticWellEq
69 unsigned int num_blocks = 0; // total number of blocks in all wells
70 unsigned int num_std_wells = 0; // number of StandardWells in this object
71 unsigned int num_ms_wells = 0; // number of MultisegmentWells in this object, must equal multisegments.size()
72 unsigned int num_blocks_so_far = 0; // keep track of where next data is written
73 unsigned int num_std_wells_so_far = 0; // keep track of where next data is written
74 std::vector<unsigned int> val_pointers; // val_pointers[wellID] == index of first block for this well in Ccols and Bcols
75
76 std::vector<std::unique_ptr<MultisegmentWellContribution<Scalar>>> multisegments;
77
78public:
79 unsigned int getNumWells(){
81 }
82
85 void addNumBlocks(unsigned int numBlocks);
86
88 void alloc();
89
92
96 void setBlockSize(unsigned int dim, unsigned int dim_wells);
97
100 void setVectorSize(unsigned N);
101
107 void addMatrix(MatrixType type, int* colIndices, Scalar* values, unsigned int val_size);
108
123 unsigned int dim_wells,
124 unsigned int Mb,
125 std::vector<Scalar>& Bvalues,
126 std::vector<unsigned int>& BcolIndices,
127 std::vector<unsigned int>& BrowPointers,
128 unsigned int DnumBlocks,
129 Scalar* Dvalues,
130 UMFPackIndex* DcolPointers,
131 UMFPackIndex* DrowIndices,
132 std::vector<Scalar>& Cvalues);
133protected:
135 virtual void APIalloc() {}
136
138 virtual void APIaddMatrix(MatrixType, int*, Scalar*, unsigned int) {}
139};
140
141} //namespace Opm
142
143#endif
Definition: WellContributions.hpp:51
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:68
void addNumBlocks(unsigned int numBlocks)
SuiteSparse_long UMFPackIndex
Definition: WellContributions.hpp:55
unsigned int num_std_wells_so_far
Definition: WellContributions.hpp:73
void setVectorSize(unsigned N)
unsigned int N
Definition: WellContributions.hpp:66
virtual void APIaddMatrix(MatrixType, int *, Scalar *, unsigned int)
Api specific upload of matrix.
Definition: WellContributions.hpp:138
unsigned int num_ms_wells
Definition: WellContributions.hpp:71
unsigned int num_std_wells
Definition: WellContributions.hpp:70
unsigned int num_blocks
Definition: WellContributions.hpp:69
bool allocated
Definition: WellContributions.hpp:64
MatrixType
StandardWell has C, D and B matrices that need to be copied.
Definition: WellContributions.hpp:57
std::vector< std::unique_ptr< MultisegmentWellContribution< Scalar > > > multisegments
Definition: WellContributions.hpp:76
unsigned int num_blocks_so_far
Definition: WellContributions.hpp:72
static std::unique_ptr< WellContributions > create(const std::string &accelerator_mode, bool useWellConn)
unsigned int getNumWells()
Definition: WellContributions.hpp:79
virtual ~WellContributions()
Empty destructor.
unsigned int dim
Definition: WellContributions.hpp:67
virtual void APIalloc()
API specific allocation.
Definition: WellContributions.hpp:135
std::vector< unsigned int > val_pointers
Definition: WellContributions.hpp:74
Definition: blackoilboundaryratevector.hh:39