[Ebos,bugfix] Correctly compute the global number of cells.

Currently, all parallel DUNE grid store some cells in addition to
interior cells. Therefore assuming that the global number of cells
(i.e. the number of cells a sequential grid needs to cover the same
whole domain with indentical cells) is not the sum of the number of
cells of the local grid. Previously, the latter was used.
This commit is contained in:
Markus Blatt 2017-03-17 10:26:00 +01:00
parent 5517b8a3f3
commit bb7934b1a2
2 changed files with 46 additions and 2 deletions

View File

@ -156,6 +156,51 @@ namespace detail {
}
}
/// \brief Get the number of local interior cells in a grid.
/// \tparam The type of the DUNE grid.
/// \param grid The grid which cells we count
/// \return The number of interior cell in the partition of the
/// grid stored on this process.
template<class Grid>
std::size_t countLocalInteriorCells(const Grid& grid)
{
if ( grid.comm().size() == 1)
{
return grid.size(0);
}
std::size_t count = 0;
const auto& gridView = grid.leafGridView();
for(auto cell = gridView.template begin<0>(),
endCell = gridView.template end<0>();
cell != endCell; ++cell)
{
if ( cell.partitionType() == Dune::InteriorEntity )
{
++count;
}
}
return count;
}
/// \brief Get the number of cells of a global grid.
///
/// In a parallel run this is the number of cells that a grid would
/// have if the whole grid was stored on one process only.
/// \tparam The type of the DUNE grid.
/// \param grid The grid which cells we count
/// \return The global number of cells.
template<class Grid>
std::size_t countGlobalCells(const Grid& grid)
{
if ( grid.comm().size() == 1)
{
return grid.size(0);
}
std::size_t count = countLocalInteriorCells(grid);
return grid.comm().sum(count);
}
template <class Scalar>
inline
double

View File

@ -190,9 +190,8 @@ namespace Opm {
const std::vector<double> depth(geo_.z().data(), geo_.z().data() + geo_.z().size());
well_model_.init(fluid_.phaseUsage(), active_, &vfp_properties_, gravity, depth, pv, &rate_converter_);
wellModel().setWellsActive( localWellsActive() );
global_nc_ = Opm::AutoDiffGrid::numCells(grid_);
// compute global sum of number of cells
global_nc_ = grid_.comm().sum( global_nc_ );
global_nc_ = detail::countGlobalCells(grid_);
if (!istlSolver_)
{