MultisegmentWellContribution.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 MULTISEGMENTWELLCONTRIBUTION_HEADER_INCLUDED
21#define MULTISEGMENTWELLCONTRIBUTION_HEADER_INCLUDED
22
23#include <vector>
24
25#if HAVE_CUDA
26#include <cuda_runtime.h>
27#endif
28
29#include <umfpack.h>
30#include <dune/common/version.hh>
31
32namespace Opm
33{
34
41
42template<class Scalar>
44{
45
46private:
47 unsigned int dim; // size of blockvectors in vectors x and y, equal to MultisegmentWell::numEq
48 unsigned int dim_wells; // size of blocks in C, B and D, equal to MultisegmentWell::numWellEq
49 unsigned int M; // number of rows, M == dim_wells*Mb
50 unsigned int Mb; // number of blockrows in C, D and B
51
52#if HAVE_CUDA
53 cudaStream_t stream; // not actually used yet, will be when MultisegmentWellContribution are applied on GPU
54#endif
55
56 // C and B are stored in BCRS format, D is stored in CSC format (Dune::UMFPack)
57 // Sparsity pattern for C is not stored, since it is the same as B
58 unsigned int DnumBlocks; // number of blocks in D
59 std::vector<Scalar> Cvals;
60 std::vector<Scalar> Dvals;
61 std::vector<Scalar> Bvals;
62 std::vector<int> Dcols; // Columnpointers, contains M+1 entries
63 std::vector<unsigned int> Bcols;
64 std::vector<int> Drows; // Rowindicies, contains DnumBlocks*dim*dim_wells entries
65 std::vector<unsigned int> Brows;
66 std::vector<Scalar> z1; // z1 = B * x
67 std::vector<Scalar> z2; // z2 = D^-1 * B * x
68 void *UMFPACK_Symbolic, *UMFPACK_Numeric;
69
72 unsigned int getColIdx(unsigned int idx);
73
74public:
75 using UMFPackIndex = SuiteSparse_long;
76
77#if HAVE_CUDA
80 void setCudaStream(cudaStream_t stream);
81#endif
82
97 MultisegmentWellContribution(unsigned int dim, unsigned int dim_wells,
98 unsigned int Mb,
99 std::vector<Scalar>& Bvalues,
100 std::vector<unsigned int>& BcolIndices,
101 std::vector<unsigned int>& BrowPointers,
102 unsigned int DnumBlocks,
103 Scalar* Dvalues,
104 UMFPackIndex* DcolPointers,
105 UMFPackIndex* DrowIndices,
106 std::vector<Scalar>& Cvalues);
107
110
115 void apply(Scalar* h_x, Scalar* h_y);
116};
117
118} //namespace Opm
119
120#endif
Definition: MultisegmentWellContribution.hpp:44
SuiteSparse_long UMFPackIndex
Definition: MultisegmentWellContribution.hpp:75
MultisegmentWellContribution(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)
~MultisegmentWellContribution()
Destroy a MultisegmentWellContribution, and free memory.
void apply(Scalar *h_x, Scalar *h_y)
Definition: blackoilboundaryratevector.hh:39