2016-06-07 05:57:25 -05:00
|
|
|
/*
|
|
|
|
Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
|
|
|
|
Copyright 2014, 2015 Dr. Blatt - HPC-Simulation-Software & Services
|
|
|
|
Copyright 2014, 2015 Statoil ASA.
|
|
|
|
Copyright 2015 NTNU
|
|
|
|
Copyright 2015 IRIS AS
|
|
|
|
|
|
|
|
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/>.
|
|
|
|
*/
|
|
|
|
|
2021-04-29 02:45:06 -05:00
|
|
|
#ifndef OPM_COUNTGLOBALCELLS_HEADER_INCLUDED
|
|
|
|
#define OPM_COUNTGLOBALCELLS_HEADER_INCLUDED
|
2016-06-07 05:57:25 -05:00
|
|
|
|
2024-06-25 03:20:41 -05:00
|
|
|
#include <opm/simulators/utils/BlackoilPhases.hpp>
|
2016-06-07 05:57:25 -05:00
|
|
|
|
2021-04-29 02:45:06 -05:00
|
|
|
#include <dune/grid/common/gridview.hh>
|
|
|
|
|
2023-09-12 05:37:10 -05:00
|
|
|
#include <algorithm>
|
2021-04-29 02:45:06 -05:00
|
|
|
#include <vector>
|
2020-02-19 03:14:49 -06:00
|
|
|
|
2016-06-07 05:57:25 -05:00
|
|
|
namespace Opm {
|
|
|
|
namespace detail {
|
|
|
|
|
2023-09-12 05:37:10 -05:00
|
|
|
/// \brief Get the number of local interior cells in a grid view.
|
|
|
|
///
|
|
|
|
/// \tparam GridView Type of DUNE grid view.
|
|
|
|
///
|
|
|
|
/// \param[in] gridView Grid view for which to count the number of
|
|
|
|
/// interior cells.
|
|
|
|
///
|
|
|
|
/// \return The number of interior cell in the partition of the grid
|
|
|
|
/// stored on this process.
|
|
|
|
template <class GridView>
|
|
|
|
std::size_t countLocalInteriorCellsGridView(const GridView& gridView)
|
|
|
|
{
|
|
|
|
if (gridView.comm().size() == 1) {
|
|
|
|
return gridView.size(0);
|
|
|
|
}
|
|
|
|
|
|
|
|
return std::distance(gridView.template begin<0, Dune::Interior_Partition>(),
|
|
|
|
gridView.template end<0, Dune::Interior_Partition>());
|
|
|
|
}
|
|
|
|
|
2017-03-17 04:26:00 -05:00
|
|
|
/// \brief Get the number of local interior cells in a grid.
|
2023-09-12 05:37:10 -05:00
|
|
|
///
|
|
|
|
/// \tparam Grid Type of DUNE grid.
|
|
|
|
///
|
|
|
|
/// \param[in] grid Grid for which to count interior cells.
|
|
|
|
///
|
|
|
|
/// \return The number of interior cells in the partition of the
|
|
|
|
/// grid stored on this process.
|
2017-03-17 04:26:00 -05:00
|
|
|
template<class Grid>
|
|
|
|
std::size_t countLocalInteriorCells(const Grid& grid)
|
|
|
|
{
|
2023-09-12 05:37:10 -05:00
|
|
|
return countLocalInteriorCellsGridView(grid.leafGridView());
|
2017-03-17 04:26:00 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
/// \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.
|
2023-09-12 05:37:10 -05:00
|
|
|
///
|
|
|
|
/// \tparam Grid Type of DUNE grid.
|
|
|
|
///
|
|
|
|
/// \param[in] grid Grid for which to count global cells.
|
|
|
|
///
|
2017-03-17 04:26:00 -05:00
|
|
|
/// \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);
|
|
|
|
}
|
2023-09-12 05:37:10 -05:00
|
|
|
const std::size_t count = countLocalInteriorCellsGridView(grid.leafGridView());
|
2017-03-17 04:26:00 -05:00
|
|
|
return grid.comm().sum(count);
|
|
|
|
}
|
|
|
|
|
2016-06-07 05:57:25 -05:00
|
|
|
} // namespace detail
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
#endif // OPM_BLACKOILDETAILS_HEADER_INCLUDED
|