My Project
WellContributions.hpp
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 <opm/simulators/linalg/bda/MultisegmentWellContribution.hpp>
27 #if HAVE_SUITESPARSE_UMFPACK
28 #include<umfpack.h>
29 #endif
30 #include <dune/common/version.hh>
31 
32 namespace Opm
33 {
34 
41 
53 {
54 public:
55  static std::unique_ptr<WellContributions> create(const std::string& accelerator_mode, bool useWellConn);
56 
57 #if DUNE_VERSION_NEWER(DUNE_ISTL, 2, 7)
58  using UMFPackIndex = SuiteSparse_long;
59 #else
60  using UMFPackIndex = int;
61 #endif
63  enum class MatrixType {
64  C,
65  D,
66  B
67  };
68 
69 protected:
70  bool allocated = false;
71 
72  unsigned int N; // number of rows (not blockrows) in vectors x and y
73  unsigned int dim; // number of columns in blocks in B and C, equal to StandardWell::numEq
74  unsigned int dim_wells; // number of rows in blocks in B and C, equal to StandardWell::numStaticWellEq
75  unsigned int num_blocks = 0; // total number of blocks in all wells
76  unsigned int num_std_wells = 0; // number of StandardWells in this object
77  unsigned int num_ms_wells = 0; // number of MultisegmentWells in this object, must equal multisegments.size()
78  unsigned int num_blocks_so_far = 0; // keep track of where next data is written
79  unsigned int num_std_wells_so_far = 0; // keep track of where next data is written
80  std::vector<unsigned int> val_pointers; // val_pointers[wellID] == index of first block for this well in Ccols and Bcols
81 
82  std::vector<std::unique_ptr<MultisegmentWellContribution>> multisegments;
83 
84 public:
85  unsigned int getNumWells(){
86  return num_std_wells + num_ms_wells;
87  }
88 
91  void addNumBlocks(unsigned int numBlocks);
92 
94  void alloc();
95 
97  virtual ~WellContributions() = default;
98 
102  void setBlockSize(unsigned int dim, unsigned int dim_wells);
103 
109  void addMatrix(MatrixType type, int *colIndices, double *values, unsigned int val_size);
110 
124  void addMultisegmentWellContribution(unsigned int dim, unsigned int dim_wells,
125  unsigned int Mb,
126  std::vector<double> &Bvalues, std::vector<unsigned int> &BcolIndices, std::vector<unsigned int> &BrowPointers,
127  unsigned int DnumBlocks, double *Dvalues,
128  UMFPackIndex *DcolPointers, UMFPackIndex *DrowIndices,
129  std::vector<double> &Cvalues);
130 protected:
132  virtual void APIalloc() {}
133 
135  virtual void APIaddMatrix(MatrixType, int*, double*, unsigned int) {}
136 };
137 } //namespace Opm
138 
139 #endif
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition: WellContributions.hpp:53
virtual ~WellContributions()=default
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:92
virtual void APIalloc()
API specific allocation.
Definition: WellContributions.hpp:132
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:71
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:123
MatrixType
StandardWell has C, D and B matrices that need to be copied.
Definition: WellContributions.hpp:63
void alloc()
Allocate memory for the StandardWells.
Definition: WellContributions.cpp:113
void addNumBlocks(unsigned int numBlocks)
Indicate how large the next StandardWell is, this function cannot be called after alloc() is called.
Definition: WellContributions.cpp:104
virtual void APIaddMatrix(MatrixType, int *, double *, unsigned int)
Api specific upload of matrix.
Definition: WellContributions.hpp:135
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27