11 #ifndef vtk_m_filter_flow_internal_BoundsMap_h 
   12 #define vtk_m_filter_flow_internal_BoundsMap_h 
   22 #include <vtkm/thirdparty/diy/diy.h> 
   24 #ifdef VTKM_ENABLE_MPI 
   26 #include <vtkm/thirdparty/diy/mpi-cast.h> 
   52     this->Init({ dataSet }, { blockId });
 
   55   BoundsMap(
const std::vector<vtkm::cont::DataSet>& dataSets) { this->Init(dataSets); }
 
   64   vtkm::Bounds GetGlobalBounds()
 const { 
return this->GlobalBounds; }
 
   68     VTKM_ASSERT(idx >= 0 && 
static_cast<std::size_t
>(idx) < this->BlockBounds.size());
 
   70     return this->BlockBounds[
static_cast<std::size_t
>(idx)];
 
   75     VTKM_ASSERT(idx >= 0 && idx < this->LocalNumBlocks);
 
   76     return this->LocalIDs[
static_cast<std::size_t
>(idx)];
 
   79   std::vector<int> FindRank(
vtkm::Id blockId)
 const 
   81     auto it = this->BlockToRankMap.find(blockId);
 
   82     if (it == this->BlockToRankMap.end())
 
   88   std::vector<vtkm::Id> FindBlocks(
const vtkm::Vec3f& p)
 const { 
return this->FindBlocks(p, -1); }
 
   90   std::vector<vtkm::Id> FindBlocks(
const vtkm::Vec3f& p,
 
   91                                    const std::vector<vtkm::Id>& ignoreBlocks)
 const 
   93     vtkm::Id ignoreID = (ignoreBlocks.empty() ? -1 : ignoreBlocks[0]);
 
   94     return FindBlocks(p, ignoreID);
 
   99     std::vector<vtkm::Id> blockIDs;
 
  100     if (this->GlobalBounds.Contains(p))
 
  103       for (
auto& it : this->BlockBounds)
 
  105         if (blockId != ignoreBlock && it.Contains(p))
 
  106           blockIDs.emplace_back(blockId);
 
  114   vtkm::Id GetTotalNumBlocks()
 const { 
return this->TotalNumBlocks; }
 
  115   vtkm::Id GetLocalNumBlocks()
 const { 
return this->LocalNumBlocks; }
 
  118   void Init(
const std::vector<vtkm::cont::DataSet>& dataSets, 
const std::vector<vtkm::Id>& blockIds)
 
  120     if (dataSets.size() != blockIds.size())
 
  123     this->LocalIDs = blockIds;
 
  124     this->LocalNumBlocks = dataSets.size();
 
  129     vtkm::Id locMinId = 0, locMaxId = 0;
 
  130     if (!this->LocalIDs.empty())
 
  132       locMinId = *std::min_element(this->LocalIDs.begin(), this->LocalIDs.end());
 
  133       locMaxId = *std::max_element(this->LocalIDs.begin(), this->LocalIDs.end());
 
  136     vtkm::Id globalMinId = 0, globalMaxId = 0;
 
  138     vtkmdiy::mpi::all_reduce(comm, locMinId, globalMinId, vtkmdiy::mpi::minimum<vtkm::Id>{});
 
  139     vtkmdiy::mpi::all_reduce(comm, locMaxId, globalMaxId, vtkmdiy::mpi::maximum<vtkm::Id>{});
 
  140     if (globalMinId != 0 || (globalMaxId - globalMinId) < 1)
 
  144     std::vector<vtkm::Id> locBlockCounts(comm.size(), 0), globalBlockCounts(comm.size(), 0);
 
  145     locBlockCounts[comm.rank()] = this->LocalIDs.size();
 
  146     vtkmdiy::mpi::all_reduce(comm, locBlockCounts, globalBlockCounts, std::plus<vtkm::Id>{});
 
  150       std::accumulate(globalBlockCounts.begin(), globalBlockCounts.end(), 
vtkm::Id{ 0 });
 
  154     for (
int i = 0; i < comm.rank(); i++)
 
  155       offset += globalBlockCounts[i];
 
  158     std::vector<vtkm::Id> localBlockIds(globalNumBlocks, 0);
 
  160     for (
const auto& bid : this->LocalIDs)
 
  161       localBlockIds[offset + idx++] = bid;
 
  164     std::vector<vtkm::Id> globalBlockIds(globalNumBlocks, 0);
 
  165     vtkmdiy::mpi::all_reduce(comm, localBlockIds, globalBlockIds, std::plus<vtkm::Id>{});
 
  170     std::vector<std::vector<vtkm::Id>> rankToBlockIds(comm.size());
 
  173     for (
int rank = 0; rank < comm.size(); rank++)
 
  175       vtkm::Id numBIds = globalBlockCounts[rank];
 
  176       rankToBlockIds[rank].resize(numBIds);
 
  177       for (
vtkm::Id i = 0; i < numBIds; i++)
 
  178         rankToBlockIds[rank][i] = globalBlockIds[offset + i];
 
  184     std::set<vtkm::Id> globalUniqueBlockIds;
 
  185     globalUniqueBlockIds.insert(globalBlockIds.begin(), globalBlockIds.end());
 
  186     this->TotalNumBlocks = globalUniqueBlockIds.size();
 
  189     std::vector<std::vector<vtkm::Id>> blockIdsToRank(this->TotalNumBlocks);
 
  190     for (
int rank = 0; rank < comm.size(); rank++)
 
  192       for (
const auto& bid : rankToBlockIds[rank])
 
  194         blockIdsToRank[bid].push_back(rank);
 
  195         this->BlockToRankMap[bid].push_back(rank);
 
  199     this->Build(dataSets);
 
  202   void Init(
const std::vector<vtkm::cont::DataSet>& dataSets)
 
  204     this->LocalNumBlocks = dataSets.size();
 
  207     this->TotalNumBlocks = assigner.nblocks();
 
  208     std::vector<int> ids;
 
  211     assigner.local_gids(Comm.rank(), ids);
 
  212     for (
const auto& i : ids)
 
  213       this->LocalIDs.emplace_back(
static_cast<vtkm::Id>(i));
 
  215     for (
vtkm::Id id = 0; 
id < this->TotalNumBlocks; 
id++)
 
  216       this->BlockToRankMap[
id] = { assigner.rank(
static_cast<int>(
id)) };
 
  217     this->Build(dataSets);
 
  220   void Build(
const std::vector<vtkm::cont::DataSet>& dataSets)
 
  222     std::vector<vtkm::Float64> vals(
static_cast<std::size_t
>(this->TotalNumBlocks * 6), 0);
 
  223     std::vector<vtkm::Float64> vals2(vals.size());
 
  225     std::vector<vtkm::Float64> localMins((this->TotalNumBlocks * 3),
 
  226                                          std::numeric_limits<vtkm::Float64>::max());
 
  227     std::vector<vtkm::Float64> localMaxs((this->TotalNumBlocks * 3),
 
  228                                          -std::numeric_limits<vtkm::Float64>::max());
 
  230     for (std::size_t i = 0; i < this->LocalIDs.size(); i++)
 
  235       vtkm::Id localID = this->LocalIDs[i];
 
  236       localMins[localID * 3 + 0] = bounds.
X.
Min;
 
  237       localMins[localID * 3 + 1] = bounds.Y.Min;
 
  238       localMins[localID * 3 + 2] = bounds.Z.Min;
 
  239       localMaxs[localID * 3 + 0] = bounds.X.Max;
 
  240       localMaxs[localID * 3 + 1] = bounds.Y.Max;
 
  241       localMaxs[localID * 3 + 2] = bounds.Z.Max;
 
  244     std::vector<vtkm::Float64> globalMins, globalMaxs;
 
  246 #ifdef VTKM_ENABLE_MPI 
  247     globalMins.resize(this->TotalNumBlocks * 3);
 
  248     globalMaxs.resize(this->TotalNumBlocks * 3);
 
  252     vtkmdiy::mpi::all_reduce(comm, localMins, globalMins, vtkmdiy::mpi::minimum<vtkm::Float64>{});
 
  253     vtkmdiy::mpi::all_reduce(comm, localMaxs, globalMaxs, vtkmdiy::mpi::maximum<vtkm::Float64>{});
 
  255     globalMins = localMins;
 
  256     globalMaxs = localMaxs;
 
  259     this->BlockBounds.resize(
static_cast<std::size_t
>(this->TotalNumBlocks));
 
  263     for (
auto& block : this->BlockBounds)
 
  270                            globalMaxs[idx + 2]);
 
  271       this->GlobalBounds.Include(block);
 
  277   std::vector<vtkm::Id> LocalIDs;
 
  278   std::map<vtkm::Id, std::vector<vtkm::Int32>> BlockToRankMap;
 
  280   std::vector<vtkm::Bounds> BlockBounds;
 
  289 #endif //vtk_m_filter_flow_internal_BoundsMap_h