My Project
findOverlapRowsAndColumns.hpp
1 /*
2  Copyright 2018 Andreas Thune
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 OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
21 #define OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
22 
23 #include <vector>
24 #include <utility>
25 #include <opm/grid/common/WellConnections.hpp>
26 
27 namespace Opm
28 {
29 namespace detail
30 {
41  template<class Grid, class W>
42  void setWellConnections(const Grid& grid, const W& wells, bool useWellConn, std::vector<std::set<int>>& wellGraph)
43  {
44  if ( grid.comm().size() > 1)
45  {
46  Dune::CartesianIndexMapper< Grid > cartMapper( grid );
47  const int numCells = cartMapper.compressedSize(); // grid.numCells()
48  wellGraph.resize(numCells);
49 
50  if (useWellConn) {
51  const auto& cpgdim = cartMapper.cartesianDimensions();
52 
53  std::vector<int> cart(cpgdim[0]*cpgdim[1]*cpgdim[2], -1);
54 
55  for( int i=0; i < numCells; ++i )
56  cart[ cartMapper.cartesianIndex( i ) ] = i;
57 
58  Dune::cpgrid::WellConnections well_indices;
59  well_indices.init(wells, cpgdim, cart);
60 
61  for (auto& well : well_indices)
62  {
63  for (auto perf = well.begin(); perf != well.end(); ++perf)
64  {
65  auto perf2 = perf;
66  for (++perf2; perf2 != well.end(); ++perf2)
67  {
68  wellGraph[*perf].insert(*perf2);
69  wellGraph[*perf2].insert(*perf);
70  }
71  }
72  }
73  }
74  }
75  }
76 
85  template<class Grid, class Mapper>
86  void findOverlapAndInterior(const Grid& grid, const Mapper& mapper, std::vector<int>& overlapRows,
87  std::vector<int>& interiorRows)
88  {
89  //only relevant in parallel case.
90  if ( grid.comm().size() > 1)
91  {
92  //Numbering of cells
93  const auto& gridView = grid.leafGridView();
94  auto elemIt = gridView.template begin<0>();
95  const auto& elemEndIt = gridView.template end<0>();
96 
97  //loop over cells in mesh
98  for (; elemIt != elemEndIt; ++elemIt)
99  {
100  const auto& elem = *elemIt;
101  int lcell = mapper.index(elem);
102 
103  if (elem.partitionType() != Dune::InteriorEntity)
104  {
105  //add row to list
106  overlapRows.push_back(lcell);
107  } else {
108  interiorRows.push_back(lcell);
109  }
110  }
111  }
112  }
113 
119  template <class Grid>
120  size_t numMatrixRowsToUseInSolver(const Grid& grid, bool ownerFirst)
121  {
122  size_t numInterior = 0;
123  if (!ownerFirst || grid.comm().size()==1)
124  return grid.leafGridView().size(0);
125  const auto& gridView = grid.leafGridView();
126  auto elemIt = gridView.template begin<0>();
127  const auto& elemEndIt = gridView.template end<0>();
128 
129  // loop over cells in mesh
130  for (; elemIt != elemEndIt; ++elemIt) {
131 
132  // Count only the interior cells.
133  if (elemIt->partitionType() == Dune::InteriorEntity) {
134  numInterior++;
135  }
136  }
137 
138  return numInterior;
139  }
140 } // namespace detail
141 } // namespace Opm
142 
143 #endif // OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27