findOverlapRowsAndColumns.hpp
Go to the documentation of this file.
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 <opm/grid/common/WellConnections.hpp>
24#include <opm/grid/common/CartesianIndexMapper.hpp>
25
26#include <cstddef>
27#include <utility>
28#include <vector>
29
30namespace Dune {
31template<class Grid> class CartesianIndexMapper;
32}
33
34namespace Opm
35{
36namespace detail
37{
48 template<class Grid, class CartMapper, class W>
49 void setWellConnections(const Grid& grid, const CartMapper& cartMapper, const W& wells, bool useWellConn, std::vector<std::set<int>>& wellGraph, int numJacobiBlocks)
50 {
51 if ( grid.comm().size() > 1 || numJacobiBlocks > 1)
52 {
53 const int numCells = cartMapper.compressedSize(); // grid.numCells()
54 wellGraph.resize(numCells);
55
56 if (useWellConn) {
57 const auto& cpgdim = cartMapper.cartesianDimensions();
58
59 std::vector<int> cart(cpgdim[0]*cpgdim[1]*cpgdim[2], -1);
60
61 for( int i=0; i < numCells; ++i )
62 cart[ cartMapper.cartesianIndex( i ) ] = i;
63
64 Dune::cpgrid::WellConnections well_indices;
65 well_indices.init(wells, cpgdim, cart);
66
67 for (auto& well : well_indices)
68 {
69 for (auto perf = well.begin(); perf != well.end(); ++perf)
70 {
71 auto perf2 = perf;
72 for (++perf2; perf2 != well.end(); ++perf2)
73 {
74 wellGraph[*perf].insert(*perf2);
75 wellGraph[*perf2].insert(*perf);
76 }
77 }
78 }
79 }
80 }
81 }
82
91 template<class Grid, class Mapper>
92 void findOverlapAndInterior(const Grid& grid, const Mapper& mapper, std::vector<int>& overlapRows,
93 std::vector<int>& interiorRows)
94 {
95 //only relevant in parallel case.
96 if ( grid.comm().size() > 1)
97 {
98 //Numbering of cells
99 const auto& gridView = grid.leafGridView();
100 //loop over cells in mesh
101 for (const auto& elem : elements(gridView))
102 {
103 int lcell = mapper.index(elem);
104
105 if (elem.partitionType() != Dune::InteriorEntity)
106 {
107 //add row to list
108 overlapRows.push_back(lcell);
109 } else {
110 interiorRows.push_back(lcell);
111 }
112 }
113 }
114 }
115
121 template <class Grid>
122 std::size_t numMatrixRowsToUseInSolver(const Grid& grid, bool ownerFirst)
123 {
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();
128
129 // loop over cells in mesh
130 const auto& range = elements(gridView, Dune::Partitions::interior);
131 numInterior = std::distance(range.begin(), range.end());
132
133 return numInterior;
134 }
135} // namespace detail
136} // namespace Opm
137
138#endif // OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
Definition: SupportsFaceTag.hpp:27
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, bool useWellConn, std::vector< std::set< int > > &wellGraph, int numJacobiBlocks)
Find cell IDs for wells contained in local grid.
Definition: findOverlapRowsAndColumns.hpp:49
Definition: BlackoilPhases.hpp:27