My Project
Loading...
Searching...
No Matches
GridPartitioning.hpp
1//===========================================================================
2//
3// File: GridPartitioning.hpp
4//
5// Created: Mon Sep 7 10:09:13 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// Bård Skaflestad <bard.skaflestad@sintef.no>
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17 Copyright 2009, 2010, 2013 SINTEF ICT, Applied Mathematics.
18 Copyright 2009, 2010 Statoil ASA.
19 Copyright 2013 Dr. Markus Blatt - HPC-Simulation-Software & Services
20 Copyright 2020, 2021 OPM-OP AS
21
22 This file is part of The Open Porous Media project (OPM).
23
24 OPM is free software: you can redistribute it and/or modify
25 it under the terms of the GNU General Public License as published by
26 the Free Software Foundation, either version 3 of the License, or
27 (at your option) any later version.
28
29 OPM is distributed in the hope that it will be useful,
30 but WITHOUT ANY WARRANTY; without even the implied warranty of
31 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
32 GNU General Public License for more details.
33
34 You should have received a copy of the GNU General Public License
35 along with OPM. If not, see <http://www.gnu.org/licenses/>.
36*/
37
38#ifndef OPM_GRIDPARTITIONING_HEADER
39#define OPM_GRIDPARTITIONING_HEADER
40
41#include <vector>
42#include <array>
43#include <set>
44#include <tuple>
45
46#include <dune/common/parallel/mpihelper.hh>
47
48#include <opm/grid/utility/OpmWellType.hpp>
49#include <opm/grid/common/WellConnections.hpp>
50#include <opm/grid/common/ZoltanGraphFunctions.hpp>
51
52namespace Dune
53{
54
55 class CpGrid;
56
58 {
59 bool operator()(const std::pair<int,int>& o, const std::pair<int,int>& v)
60 {
61 return o.first < v.first;
62 }
63 };
64
73 void partition(const CpGrid& grid,
74 const std::array<int, 3>& initial_split,
75 int& num_part,
76 std::vector<int>& cell_part,
77 bool recursive = false,
78 bool ensureConnectivity = true);
79
88 void addOverlapLayer(const CpGrid& grid,
89 const std::vector<int>& cell_part,
90 std::vector<std::set<int> >& cell_overlap,
91 int mypart, int overlapLayers, bool all=false);
92
105 int addOverlapLayer(const CpGrid& grid, const std::vector<int>& cell_part,
106 std::vector<std::tuple<int,int,char>>& exportList,
107 std::vector<std::tuple<int,int,char,int>>& importList,
108 const Communication<Dune::MPIHelper::MPICommunicator>& cc,
109 bool addCornerCells, const double* trans, int layers = 1);
110
111namespace cpgrid
112{
113#if HAVE_MPI
135 std::tuple<std::vector<int>, std::vector<std::pair<std::string,bool>>,
136 std::vector<std::tuple<int,int,char> >,
137 std::vector<std::tuple<int,int,char,int> >,
138 WellConnections>
139 createListsFromParts(const CpGrid& grid, const std::vector<cpgrid::OpmWellType> * wells,
140 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
141 const double* transmissibilities, const std::vector<int>& parts,
142 bool allowDistributedWells, std::shared_ptr<cpgrid::CombinedGridWellGraph> gridAndWells = nullptr);
143
164 std::tuple<std::vector<int>, std::vector<std::pair<std::string,bool>>,
165 std::vector<std::tuple<int,int,char> >,
166 std::vector<std::tuple<int,int,char,int> >,
167 WellConnections>
168 vanillaPartitionGridOnRoot(const CpGrid& grid,
169 const std::vector<cpgrid::OpmWellType> * wells,
170 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
171 const double* transmissibilities,
172 bool allowDistributedWells);
173#endif
174} // namespace cpgrid
175} // namespace Dune
176
177
178#endif // OPM_GRIDPARTITIONING_HEADER
[ provides Dune::Grid ]
Definition CpGrid.hpp:201
Copyright 2019 Equinor AS.
Definition CartesianIndexMapper.hpp:10
void partition(const CpGrid &grid, const coord_t &initial_split, int &num_part, std::vector< int > &cell_part, bool recursive, bool ensureConnectivity)
Partition a CpGrid based on (ijk) coordinates, with splitting to ensure that each partition is connec...
Definition GridPartitioning.cpp:195
Definition GridPartitioning.hpp:58