2019-06-19 08:04:14 -05:00
|
|
|
/*
|
|
|
|
Copyright 2018 Andreas Thune
|
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
*/
|
|
|
|
|
|
|
|
#ifndef OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
|
|
|
|
#define OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
|
|
|
|
|
|
|
|
#include <vector>
|
|
|
|
#include <utility>
|
2019-11-21 10:04:45 -06:00
|
|
|
#include <opm/grid/common/WellConnections.hpp>
|
2021-12-01 07:00:21 -06:00
|
|
|
#include <opm/grid/common/CartesianIndexMapper.hpp>
|
|
|
|
|
|
|
|
namespace Dune {
|
|
|
|
template<class Grid> class CartesianIndexMapper;
|
|
|
|
}
|
2019-06-19 08:04:14 -05:00
|
|
|
|
|
|
|
namespace Opm
|
|
|
|
{
|
|
|
|
namespace detail
|
|
|
|
{
|
2019-11-21 10:04:45 -06:00
|
|
|
/// \brief Find cell IDs for wells contained in local grid.
|
|
|
|
///
|
2019-12-11 15:04:40 -06:00
|
|
|
/// Cell IDs of wells stored in a graph, so it can be used to create
|
2019-11-21 10:04:45 -06:00
|
|
|
/// an adjacency pattern. Only relevant when the UseWellContribusion option is set to true
|
|
|
|
/// \tparam The type of the DUNE grid.
|
|
|
|
/// \tparam Well vector type
|
|
|
|
/// \param grid The grid where we look for overlap cells.
|
|
|
|
/// \param wells List of wells contained in grid.
|
|
|
|
/// \param useWellConn Boolean that is true when UseWellContribusion is true
|
|
|
|
/// \param wellGraph Cell IDs of well cells stored in a graph.
|
2021-12-01 07:00:21 -06:00
|
|
|
template<class Grid, class CartMapper, class W>
|
|
|
|
void setWellConnections(const Grid& grid, const CartMapper& cartMapper, const W& wells, bool useWellConn, std::vector<std::set<int>>& wellGraph, int numJacobiBlocks)
|
2019-11-21 10:04:45 -06:00
|
|
|
{
|
2022-04-21 10:18:32 -05:00
|
|
|
if ( grid.comm().size() > 1 || numJacobiBlocks > 1)
|
2019-11-21 10:04:45 -06:00
|
|
|
{
|
2019-12-11 15:04:40 -06:00
|
|
|
const int numCells = cartMapper.compressedSize(); // grid.numCells()
|
|
|
|
wellGraph.resize(numCells);
|
2019-11-21 10:04:45 -06:00
|
|
|
|
|
|
|
if (useWellConn) {
|
2019-12-11 15:04:40 -06:00
|
|
|
const auto& cpgdim = cartMapper.cartesianDimensions();
|
2019-11-21 10:04:45 -06:00
|
|
|
|
|
|
|
std::vector<int> cart(cpgdim[0]*cpgdim[1]*cpgdim[2], -1);
|
|
|
|
|
2019-12-11 15:04:40 -06:00
|
|
|
for( int i=0; i < numCells; ++i )
|
|
|
|
cart[ cartMapper.cartesianIndex( i ) ] = i;
|
2019-11-21 10:04:45 -06:00
|
|
|
|
|
|
|
Dune::cpgrid::WellConnections well_indices;
|
|
|
|
well_indices.init(wells, cpgdim, cart);
|
|
|
|
|
|
|
|
for (auto& well : well_indices)
|
|
|
|
{
|
|
|
|
for (auto perf = well.begin(); perf != well.end(); ++perf)
|
|
|
|
{
|
|
|
|
auto perf2 = perf;
|
|
|
|
for (++perf2; perf2 != well.end(); ++perf2)
|
|
|
|
{
|
|
|
|
wellGraph[*perf].insert(*perf2);
|
|
|
|
wellGraph[*perf2].insert(*perf);
|
|
|
|
}
|
2019-12-11 15:04:40 -06:00
|
|
|
}
|
2019-11-21 10:04:45 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2019-06-19 08:04:14 -05:00
|
|
|
|
|
|
|
/// \brief Find the rows corresponding to overlap cells
|
|
|
|
///
|
2019-11-21 10:04:45 -06:00
|
|
|
/// Loop over grid and store cell ids of rows
|
2019-06-19 08:04:14 -05:00
|
|
|
/// corresponding to overlap cells.
|
|
|
|
/// \tparam The type of the DUNE grid.
|
|
|
|
/// \param grid The grid where we look for overlap cells.
|
2019-11-21 10:04:45 -06:00
|
|
|
/// \param overlapRows List where overlap rows are stored.
|
|
|
|
/// \param interiorRows List where overlap rows are stored.
|
2020-03-12 07:59:46 -05:00
|
|
|
template<class Grid, class Mapper>
|
|
|
|
void findOverlapAndInterior(const Grid& grid, const Mapper& mapper, std::vector<int>& overlapRows,
|
2019-11-21 10:04:45 -06:00
|
|
|
std::vector<int>& interiorRows)
|
2019-06-19 08:04:14 -05:00
|
|
|
{
|
2019-11-21 10:04:45 -06:00
|
|
|
//only relevant in parallel case.
|
2019-12-11 15:04:40 -06:00
|
|
|
if ( grid.comm().size() > 1)
|
2019-11-21 10:04:45 -06:00
|
|
|
{
|
|
|
|
//Numbering of cells
|
2019-06-19 08:04:14 -05:00
|
|
|
const auto& gridView = grid.leafGridView();
|
|
|
|
auto elemIt = gridView.template begin<0>();
|
|
|
|
const auto& elemEndIt = gridView.template end<0>();
|
|
|
|
|
2019-11-21 10:04:45 -06:00
|
|
|
//loop over cells in mesh
|
2019-12-11 15:04:40 -06:00
|
|
|
for (; elemIt != elemEndIt; ++elemIt)
|
|
|
|
{
|
2019-06-19 08:04:14 -05:00
|
|
|
const auto& elem = *elemIt;
|
2020-03-12 07:59:46 -05:00
|
|
|
int lcell = mapper.index(elem);
|
2019-12-11 15:04:40 -06:00
|
|
|
|
2019-11-21 10:04:45 -06:00
|
|
|
if (elem.partitionType() != Dune::InteriorEntity)
|
2019-12-11 15:04:40 -06:00
|
|
|
{
|
2019-11-21 10:04:45 -06:00
|
|
|
//add row to list
|
|
|
|
overlapRows.push_back(lcell);
|
|
|
|
} else {
|
|
|
|
interiorRows.push_back(lcell);
|
2019-06-19 08:04:14 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2020-03-08 09:45:27 -05:00
|
|
|
|
|
|
|
/// \brief If ownerFirst=true, returns the number of interior cells in grid, else just numCells().
|
|
|
|
///
|
|
|
|
/// If cells in grid is ordered so that interior/owner cells come before overlap/copy cells, the method
|
|
|
|
/// returns the number of interior cells numInterior. In the linear solver only the first numInterior rows of
|
|
|
|
/// the matrix are needed.
|
2019-12-07 08:38:19 -06:00
|
|
|
template <class Grid>
|
2020-03-08 09:45:27 -05:00
|
|
|
size_t numMatrixRowsToUseInSolver(const Grid& grid, bool ownerFirst)
|
2019-12-07 08:38:19 -06:00
|
|
|
{
|
|
|
|
size_t numInterior = 0;
|
2020-03-08 09:45:27 -05:00
|
|
|
if (!ownerFirst || grid.comm().size()==1)
|
2020-11-10 01:54:44 -06:00
|
|
|
return grid.leafGridView().size(0);
|
2019-12-07 08:38:19 -06:00
|
|
|
const auto& gridView = grid.leafGridView();
|
|
|
|
auto elemIt = gridView.template begin<0>();
|
|
|
|
const auto& elemEndIt = gridView.template end<0>();
|
|
|
|
|
|
|
|
// loop over cells in mesh
|
|
|
|
for (; elemIt != elemEndIt; ++elemIt) {
|
|
|
|
|
2020-03-08 09:45:27 -05:00
|
|
|
// Count only the interior cells.
|
|
|
|
if (elemIt->partitionType() == Dune::InteriorEntity) {
|
2019-12-07 08:38:19 -06:00
|
|
|
numInterior++;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
return numInterior;
|
|
|
|
}
|
2019-06-19 08:04:14 -05:00
|
|
|
} // namespace detail
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
#endif // OPM_FINDOVERLAPROWSANDCOLUMNS_HEADER_INCLUDED
|