//============================================================================ // Copyright (c) Kitware, Inc. // All rights reserved. // See LICENSE.txt for details. // // This software is distributed WITHOUT ANY WARRANTY; without even // the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR // PURPOSE. See the above copyright notice for more information. //============================================================================ /* * This example demonstrates how one can write a filter that uses MPI * for hybrid-parallelism. The `vtkm::filter::Histogram` is another approach for * implementing the same that uses DIY. This example doesn't use DIY, instead * uses MPI calls directly. */ #include "HistogramMPI.h" #include #include #include #include #include #include #include #include #include #include namespace { template VTKM_CONT vtkm::cont::ArrayHandle CreateArray(T min, T max, vtkm::Id numVals) { std::mt19937 gen; std::uniform_real_distribution dis(static_cast(min), static_cast(max)); vtkm::cont::ArrayHandle handle; handle.Allocate(numVals); std::generate(vtkm::cont::ArrayPortalToIteratorBegin(handle.WritePortal()), vtkm::cont::ArrayPortalToIteratorEnd(handle.WritePortal()), [&]() { return static_cast(dis(gen)); }); return handle; } } int main(int argc, char* argv[]) { //parse out all vtk-m related command line options auto opts = vtkm::cont::InitializeOptions::DefaultAnyDevice | vtkm::cont::InitializeOptions::Strict; vtkm::cont::Initialize(argc, argv, opts); // setup MPI environment. vtkmdiy::mpi::environment env(argc, argv); // will finalize on destruction vtkmdiy::mpi::communicator world; // the default is MPI_COMM_WORLD // tell VTK-m the communicator to use. vtkm::cont::EnvironmentTracker::SetCommunicator(world); int rank, size; MPI_Comm_rank(vtkmdiy::mpi::mpi_cast(world.handle()), &rank); MPI_Comm_size(vtkmdiy::mpi::mpi_cast(world.handle()), &size); if (argc != 2) { if (rank == 0) { std::cout << "Usage: " << std::endl << "$ " << argv[0] << " " << std::endl; } return EXIT_FAILURE; } const vtkm::Id num_bins = static_cast(std::atoi(argv[1])); const vtkm::Id numVals = 1024; vtkm::cont::PartitionedDataSet pds; vtkm::cont::DataSet ds; ds.AddPointField("pointvar", CreateArray(-1024, 1024, numVals)); pds.AppendPartition(ds); example::HistogramMPI histogram; histogram.SetActiveField("pointvar"); histogram.SetNumberOfBins(std::max(1, num_bins)); vtkm::cont::PartitionedDataSet result = histogram.Execute(pds); vtkm::cont::ArrayHandle bins; result.GetPartition(0).GetField("histogram").GetData().AsArrayHandle(bins); auto binPortal = bins.ReadPortal(); if (rank == 0) { // print histogram. std::cout << "Histogram (" << num_bins << ")" << std::endl; vtkm::Id count = 0; for (vtkm::Id cc = 0; cc < num_bins; ++cc) { std::cout << " bin[" << cc << "] = " << binPortal.Get(cc) << std::endl; count += binPortal.Get(cc); } if (count != numVals * size) { std::cout << "ERROR: bins mismatched!" << std::endl; return EXIT_FAILURE; } } return EXIT_SUCCESS; }