WellConnections.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 Dr. Blatt - HPC-Simulation-Software & Services.
3 Copyright 2016 Statoil AS
4
5 This file is part of The Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20#ifndef DUNE_CPGRID_WELL_CONNECTIONS_HEADER_INCLUDED
21#define DUNE_CPGRID_WELL_CONNECTIONS_HEADER_INCLUDED
22
23#include <set>
24#include <unordered_set>
25#include <vector>
26#include <functional>
27
28#ifdef HAVE_MPI
30#include "mpi.h"
32#endif
33
34#include <dune/common/parallel/communication.hh>
35
36#include <opm/grid/CpGrid.hpp>
38
39namespace Dune
40{
41namespace cpgrid
42{
43
51{
52
53public:
55 typedef std::vector<std::set<int> >::const_iterator const_iterator;
56
59
60 WellConnections() = default;
61
71 WellConnections(const std::vector<OpmWellType>& wells,
72 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
73 const std::array<int, 3>& cartesianSize,
74 const std::vector<int>& cartesian_to_compressed);
75
82 WellConnections(const std::vector<OpmWellType>& wells,
83 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
84 const Dune::CpGrid& cpGrid);
85
95 void init(const std::vector<OpmWellType>& wells,
96 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
97 const std::array<int, 3>& cartesianSize,
98 const std::vector<int>& cartesian_to_compressed);
99
104 const std::set<int>& operator[](std::size_t i) const
105 {
106 return well_indices_[i];
107 }
108
111 {
112 return well_indices_.begin();
113 }
114
117 {
118 return well_indices_.end();
119 }
120
122 std::size_t size() const
123 {
124 return well_indices_.size();
125 }
126private:
129 std::vector<std::set<int> > well_indices_;
130};
131
132
133#ifdef HAVE_MPI
148std::vector<std::vector<int> >
149perforatingWellIndicesOnProc(const std::vector<int>& parts,
150 const std::vector<Dune::cpgrid::OpmWellType>& wells,
151 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
152 const CpGrid& cpgrid);
153
171std::vector<std::vector<int> >
172postProcessPartitioningForWells(std::vector<int>& parts,
173 std::function<int(int)> gid,
174 const std::vector<OpmWellType>& wells,
175 const WellConnections& well_connections,
176 const std::vector<std::set<int> >& wellGraph,
177 std::vector<std::tuple<int,int,char>>& exportList,
178 std::vector<std::tuple<int,int,char,int>>& importList,
179 const Communication<MPI_Comm>& cc);
180
189std::vector<std::pair<std::string,bool>>
190computeParallelWells(const std::vector<std::vector<int> >& wells_on_proc,
191 const std::vector<OpmWellType>& wells,
192 const Communication<MPI_Comm>& cc,
193 int root);
194#endif
195} // end namespace cpgrid
196} // end namespace Dune
197
198#endif
[ provides Dune::Grid ]
Definition: CpGrid.hpp:198
A class calculating and representing all connections of wells.
Definition: WellConnections.hpp:51
const_iterator begin() const
Get a begin iterator.
Definition: WellConnections.hpp:110
void init(const std::vector< OpmWellType > &wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const std::array< int, 3 > &cartesianSize, const std::vector< int > &cartesian_to_compressed)
Initialze the data of the container.
const std::set< int > & operator[](std::size_t i) const
Access all connections of a well.
Definition: WellConnections.hpp:104
WellConnections(const std::vector< OpmWellType > &wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const Dune::CpGrid &cpGrid)
Constructor.
const_iterator end() const
Get the end iterator.
Definition: WellConnections.hpp:116
const_iterator iterator
The iterator type (always const).
Definition: WellConnections.hpp:58
std::vector< std::set< int > >::const_iterator const_iterator
The const iterator type.
Definition: WellConnections.hpp:55
std::size_t size() const
\breif Get the number of wells
Definition: WellConnections.hpp:122
WellConnections(const std::vector< OpmWellType > &wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const std::array< int, 3 > &cartesianSize, const std::vector< int > &cartesian_to_compressed)
Constructor.
std::vector< std::vector< int > > perforatingWellIndicesOnProc(const std::vector< int > &parts, const std::vector< Dune::cpgrid::OpmWellType > &wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const CpGrid &cpgrid)
Determines the wells that have perforate cells for each process.
std::vector< std::vector< int > > postProcessPartitioningForWells(std::vector< int > &parts, std::function< int(int)> gid, const std::vector< OpmWellType > &wells, const WellConnections &well_connections, const std::vector< std::set< int > > &wellGraph, std::vector< std::tuple< int, int, char > > &exportList, std::vector< std::tuple< int, int, char, int > > &importList, const Communication< MPI_Comm > &cc)
Computes wells assigned to processes.
std::vector< std::pair< std::string, bool > > computeParallelWells(const std::vector< std::vector< int > > &wells_on_proc, const std::vector< OpmWellType > &wells, const Communication< MPI_Comm > &cc, int root)
Computes whether wells are perforating cells on this process.
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10