mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Fixes parallel computation of average formation factors.
This was using the sum of the sum of all cells of the partitions. Hence the shared cells were counted twice.
This commit is contained in:
parent
0aa3f51eb7
commit
ca8de16285
@ -617,6 +617,31 @@ public:
|
||||
{
|
||||
return centroids_;
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Get the number of cells in the global leaf grid view.
|
||||
* \warn This is a collective operation that needs to be called
|
||||
* on all ranks.
|
||||
*/
|
||||
std::size_t noGlobalCells() const
|
||||
{
|
||||
const auto& grid = asImp_().grid();
|
||||
if (grid.comm().size() == 1)
|
||||
{
|
||||
return grid.leafGridView().size(0);
|
||||
}
|
||||
const auto& gridView = grid.leafGridView();
|
||||
const auto& elemEndIt = gridView.template end</*codim=*/0, Dune::Interior_Partition>();
|
||||
std::size_t global_nc = 0;
|
||||
|
||||
for (auto elemIt = gridView.template begin</*codim=*/0, Dune::Interior_Partition>();
|
||||
elemIt != elemEndIt; ++elemIt)
|
||||
{
|
||||
++global_nc;
|
||||
}
|
||||
return grid.comm().sum(global_nc);
|
||||
}
|
||||
|
||||
protected:
|
||||
void callImplementationInit()
|
||||
{
|
||||
|
@ -43,12 +43,9 @@ namespace Opm {
|
||||
// Create the guide rate container.
|
||||
guideRate_.reset(new GuideRate (ebosSimulator_.vanguard().schedule()));
|
||||
|
||||
// calculate the number of elements of the compressed sequential grid. this needs
|
||||
// to be done in two steps because the dune communicator expects a reference as
|
||||
// argument for sum()
|
||||
const auto& gridView = ebosSimulator_.gridView();
|
||||
number_of_cells_ = gridView.size(/*codim=*/0);
|
||||
global_nc_ = gridView.comm().sum(number_of_cells_);
|
||||
number_of_cells_ = ebosSimulator_.gridView().size(0);
|
||||
// Number of cells the global grid view
|
||||
global_nc_ = ebosSimulator_.vanguard().noGlobalCells();
|
||||
|
||||
// Set up cartesian mapping.
|
||||
const auto& grid = ebosSimulator_.vanguard().grid();
|
||||
|
Loading…
Reference in New Issue
Block a user