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 #include <opm/grid/common/CartesianIndexMapper.hpp>
27 
28 namespace Dune {
29 template<class Grid> class CartesianIndexMapper;
30 }
31 
32 namespace Opm
33 {
34 namespace detail
35 {
46  template<class Grid, class CartMapper, class W>
47  void setWellConnections(const Grid& grid, const CartMapper& cartMapper, const W& wells, bool useWellConn, std::vector<std::set<int>>& wellGraph, int numJacobiBlocks)
48  {
49  if ( grid.comm().size() > 1 || numJacobiBlocks > 1)
50  {
51  const int numCells = cartMapper.compressedSize(); // grid.numCells()
52  wellGraph.resize(numCells);
53 
54  if (useWellConn) {
55  const auto& cpgdim = cartMapper.cartesianDimensions();
56 
57  std::vector<int> cart(cpgdim[0]*cpgdim[1]*cpgdim[2], -1);
58 
59  for( int i=0; i < numCells; ++i )
60  cart[ cartMapper.cartesianIndex( i ) ] = i;
61 
62  Dune::cpgrid::WellConnections well_indices;
63  well_indices.init(wells, cpgdim, cart);
64 
65  for (auto& well : well_indices)
66  {
67  for (auto perf = well.begin(); perf != well.end(); ++perf)
68  {
69  auto perf2 = perf;
70  for (++perf2; perf2 != well.end(); ++perf2)
71  {
72  wellGraph[*perf].insert(*perf2);
73  wellGraph[*perf2].insert(*perf);
74  }
75  }
76  }
77  }
78  }
79  }
80 
89  template<class Grid, class Mapper>
90  void findOverlapAndInterior(const Grid& grid, const Mapper& mapper, std::vector<int>& overlapRows,
91  std::vector<int>& interiorRows)
92  {
93  //only relevant in parallel case.
94  if ( grid.comm().size() > 1)
95  {
96  //Numbering of cells
97  const auto& gridView = grid.leafGridView();
98  //loop over cells in mesh
99  for (const auto& elem : elements(gridView))
100  {
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 
127  // loop over cells in mesh
128  const auto& range = elements(gridView, Dune::Partitions::interior);
129  numInterior = std::distance(range.begin(), range.end());
130 
131  return numInterior;
132  }
133 } // namespace detail
134 } // namespace Opm
135 
136 #endif // OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
Definition: findOverlapRowsAndColumns.hpp:29
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27