20#ifndef OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
21#define OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
23#include <opm/grid/common/WellConnections.hpp>
24#include <opm/grid/common/CartesianIndexMapper.hpp>
31template<
class Gr
id>
class CartesianIndexMapper;
48 template<
class Gr
id,
class CartMapper,
class W>
49 void setWellConnections(
const Grid& grid,
const CartMapper& cartMapper,
const W& wells,
const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
bool useWellConn, std::vector<std::set<int>>& wellGraph,
int numJacobiBlocks)
51 if ( grid.comm().size() > 1 || numJacobiBlocks > 1)
53 const int numCells = cartMapper.compressedSize();
54 wellGraph.resize(numCells);
57 const auto& cpgdim = cartMapper.cartesianDimensions();
59 std::vector<int> cart(cpgdim[0]*cpgdim[1]*cpgdim[2], -1);
61 for(
int i=0; i < numCells; ++i )
62 cart[ cartMapper.cartesianIndex( i ) ] = i;
64 Dune::cpgrid::WellConnections well_indices;
65 well_indices.init(wells, possibleFutureConnections, cpgdim, cart);
67 for (
auto& well : well_indices)
69 for (
auto perf = well.begin(); perf != well.end(); ++perf)
72 for (++perf2; perf2 != well.end(); ++perf2)
74 wellGraph[*perf].insert(*perf2);
75 wellGraph[*perf2].insert(*perf);
91 template<
class Gr
id,
class Mapper>
93 std::vector<int>& interiorRows)
96 if ( grid.comm().size() > 1)
99 const auto& gridView = grid.leafGridView();
101 for (
const auto& elem : elements(gridView))
103 int lcell = mapper.index(elem);
105 if (elem.partitionType() != Dune::InteriorEntity)
108 overlapRows.push_back(lcell);
110 interiorRows.push_back(lcell);
121 template <
class Gr
id>
124 std::size_t numInterior = 0;
125 if (!ownerFirst || grid.comm().size()==1)
126 return grid.leafGridView().size(0);
127 const auto& gridView = grid.leafGridView();
130 const auto& range = elements(gridView, Dune::Partitions::interior);
131 numInterior = std::distance(range.begin(), range.end());
Definition: fvbaseprimaryvariables.hh:141
std::size_t numMatrixRowsToUseInSolver(const Grid &grid, bool ownerFirst)
If ownerFirst=true, returns the number of interior cells in grid, else just numCells().
Definition: findOverlapRowsAndColumns.hpp:122
void findOverlapAndInterior(const Grid &grid, const Mapper &mapper, std::vector< int > &overlapRows, std::vector< int > &interiorRows)
Find the rows corresponding to overlap cells.
Definition: findOverlapRowsAndColumns.hpp:92
void setWellConnections(const Grid &grid, const CartMapper &cartMapper, const W &wells, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, bool useWellConn, std::vector< std::set< int > > &wellGraph, int numJacobiBlocks)
Find cell IDs for wells contained in local grid.
Definition: findOverlapRowsAndColumns.hpp:49
Definition: blackoilboundaryratevector.hh:37