2018-04-11 00:11:32 +00:00
|
|
|
//============================================================================
|
|
|
|
// Copyright (c) Kitware, Inc.
|
|
|
|
// All rights reserved.
|
|
|
|
// See LICENSE.txt for details.
|
2019-04-15 23:24:21 +00:00
|
|
|
//
|
2018-04-11 00:11:32 +00:00
|
|
|
// 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.
|
|
|
|
//============================================================================
|
|
|
|
|
2022-09-26 16:23:26 +00:00
|
|
|
#include "HistogramMPI.h"
|
|
|
|
|
2022-01-17 16:18:42 +00:00
|
|
|
#include <vtkm/filter/density_estimate/worklet/FieldHistogram.h>
|
2018-04-11 00:11:32 +00:00
|
|
|
|
|
|
|
#include <vtkm/cont/Algorithm.h>
|
|
|
|
#include <vtkm/cont/ArrayPortalToIterators.h>
|
2019-09-02 14:38:47 +00:00
|
|
|
#include <vtkm/cont/AssignerPartitionedDataSet.h>
|
2018-04-11 00:11:32 +00:00
|
|
|
#include <vtkm/cont/EnvironmentTracker.h>
|
|
|
|
#include <vtkm/cont/ErrorFilterExecution.h>
|
|
|
|
#include <vtkm/cont/FieldRangeGlobalCompute.h>
|
|
|
|
|
2019-02-06 16:16:10 +00:00
|
|
|
#include <vtkm/thirdparty/diy/diy.h>
|
2020-06-08 20:57:51 +00:00
|
|
|
#include <vtkm/thirdparty/diy/mpi-cast.h>
|
|
|
|
|
|
|
|
#include <mpi.h>
|
2018-04-11 00:11:32 +00:00
|
|
|
|
|
|
|
namespace example
|
|
|
|
{
|
|
|
|
namespace detail
|
|
|
|
{
|
|
|
|
|
|
|
|
class DistributedHistogram
|
|
|
|
{
|
|
|
|
std::vector<vtkm::cont::ArrayHandle<vtkm::Id>> LocalBlocks;
|
|
|
|
|
|
|
|
public:
|
|
|
|
DistributedHistogram(vtkm::Id numLocalBlocks)
|
|
|
|
: LocalBlocks(static_cast<size_t>(numLocalBlocks))
|
|
|
|
{
|
|
|
|
}
|
|
|
|
|
|
|
|
void SetLocalHistogram(vtkm::Id index, const vtkm::cont::ArrayHandle<vtkm::Id>& bins)
|
|
|
|
{
|
|
|
|
this->LocalBlocks[static_cast<size_t>(index)] = bins;
|
|
|
|
}
|
|
|
|
|
|
|
|
void SetLocalHistogram(vtkm::Id index, const vtkm::cont::Field& field)
|
|
|
|
{
|
2021-08-06 12:53:29 +00:00
|
|
|
this->SetLocalHistogram(index,
|
|
|
|
field.GetData().AsArrayHandle<vtkm::cont::ArrayHandle<vtkm::Id>>());
|
2018-04-11 00:11:32 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id> ReduceAll(const vtkm::Id numBins) const
|
|
|
|
{
|
|
|
|
const vtkm::Id numLocalBlocks = static_cast<vtkm::Id>(this->LocalBlocks.size());
|
|
|
|
auto comm = vtkm::cont::EnvironmentTracker::GetCommunicator();
|
|
|
|
if (comm.size() == 1 && numLocalBlocks <= 1)
|
|
|
|
{
|
|
|
|
// no reduction necessary.
|
|
|
|
return numLocalBlocks == 0 ? vtkm::cont::ArrayHandle<vtkm::Id>() : this->LocalBlocks[0];
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// reduce local bins first.
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id> local;
|
|
|
|
local.Allocate(numBins);
|
2020-01-28 19:14:32 +00:00
|
|
|
std::fill(vtkm::cont::ArrayPortalToIteratorBegin(local.WritePortal()),
|
|
|
|
vtkm::cont::ArrayPortalToIteratorEnd(local.WritePortal()),
|
2018-04-11 00:11:32 +00:00
|
|
|
static_cast<vtkm::Id>(0));
|
|
|
|
for (const auto& lbins : this->LocalBlocks)
|
|
|
|
{
|
|
|
|
vtkm::cont::Algorithm::Transform(local, lbins, local, vtkm::Add());
|
|
|
|
}
|
|
|
|
|
|
|
|
// now reduce across ranks using MPI.
|
|
|
|
|
|
|
|
// converting to std::vector
|
2018-05-09 21:49:45 +00:00
|
|
|
std::vector<vtkm::Id> send_buf(static_cast<std::size_t>(numBins));
|
2020-01-28 19:14:32 +00:00
|
|
|
std::copy(vtkm::cont::ArrayPortalToIteratorBegin(local.ReadPortal()),
|
|
|
|
vtkm::cont::ArrayPortalToIteratorEnd(local.ReadPortal()),
|
2018-04-11 00:11:32 +00:00
|
|
|
send_buf.begin());
|
|
|
|
|
2018-05-09 21:49:45 +00:00
|
|
|
std::vector<vtkm::Id> recv_buf(static_cast<std::size_t>(numBins));
|
2018-04-11 00:11:32 +00:00
|
|
|
MPI_Reduce(&send_buf[0],
|
|
|
|
&recv_buf[0],
|
|
|
|
static_cast<int>(numBins),
|
|
|
|
sizeof(vtkm::Id) == 4 ? MPI_INT : MPI_LONG,
|
|
|
|
MPI_SUM,
|
|
|
|
0,
|
2020-06-08 20:57:51 +00:00
|
|
|
vtkmdiy::mpi::mpi_cast(comm.handle()));
|
2018-04-11 00:11:32 +00:00
|
|
|
|
|
|
|
if (comm.rank() == 0)
|
|
|
|
{
|
|
|
|
local.Allocate(numBins);
|
|
|
|
std::copy(recv_buf.begin(),
|
|
|
|
recv_buf.end(),
|
2020-01-28 19:14:32 +00:00
|
|
|
vtkm::cont::ArrayPortalToIteratorBegin(local.WritePortal()));
|
2018-04-11 00:11:32 +00:00
|
|
|
return local;
|
|
|
|
}
|
|
|
|
return vtkm::cont::ArrayHandle<vtkm::Id>();
|
|
|
|
}
|
|
|
|
};
|
|
|
|
|
|
|
|
} // namespace detail
|
|
|
|
|
|
|
|
//-----------------------------------------------------------------------------
|
2022-09-26 16:23:26 +00:00
|
|
|
VTKM_CONT vtkm::cont::DataSet HistogramMPI::DoExecute(const vtkm::cont::DataSet& input)
|
2018-04-11 00:11:32 +00:00
|
|
|
{
|
2022-09-26 16:23:26 +00:00
|
|
|
const auto& fieldArray = this->GetFieldFromDataSet(input).GetData();
|
|
|
|
|
|
|
|
if (!this->InExecutePartitions)
|
|
|
|
{
|
|
|
|
// Handle initialization that would be done in PreExecute if the data set had partitions.
|
|
|
|
if (this->Range.IsNonEmpty())
|
|
|
|
{
|
|
|
|
this->ComputedRange = this->Range;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
auto handle = vtkm::cont::FieldRangeGlobalCompute(
|
|
|
|
input, this->GetActiveFieldName(), this->GetActiveFieldAssociation());
|
|
|
|
if (handle.GetNumberOfValues() != 1)
|
|
|
|
{
|
|
|
|
throw vtkm::cont::ErrorFilterExecution("expecting scalar field.");
|
|
|
|
}
|
|
|
|
this->ComputedRange = handle.ReadPortal().Get(0);
|
|
|
|
}
|
|
|
|
}
|
2018-04-11 00:11:32 +00:00
|
|
|
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id> binArray;
|
|
|
|
|
2022-09-26 16:23:26 +00:00
|
|
|
auto resolveType = [&](const auto& concrete) {
|
|
|
|
using T = typename std::decay_t<decltype(concrete)>::ValueType;
|
|
|
|
T delta;
|
|
|
|
|
|
|
|
vtkm::worklet::FieldHistogram worklet;
|
|
|
|
worklet.Run(concrete,
|
2018-04-11 00:11:32 +00:00
|
|
|
this->NumberOfBins,
|
|
|
|
static_cast<T>(this->ComputedRange.Min),
|
|
|
|
static_cast<T>(this->ComputedRange.Max),
|
|
|
|
delta,
|
2018-11-29 22:02:00 +00:00
|
|
|
binArray);
|
2018-04-11 00:11:32 +00:00
|
|
|
|
2022-09-26 16:23:26 +00:00
|
|
|
this->BinDelta = static_cast<vtkm::Float64>(delta);
|
|
|
|
};
|
|
|
|
|
|
|
|
fieldArray
|
|
|
|
.CastAndCallForTypesWithFloatFallback<vtkm::TypeListFieldScalar, VTKM_DEFAULT_STORAGE_LIST>(
|
|
|
|
resolveType);
|
|
|
|
|
2018-04-11 00:11:32 +00:00
|
|
|
vtkm::cont::DataSet output;
|
2022-09-26 16:23:26 +00:00
|
|
|
output.AddField(
|
|
|
|
{ this->GetOutputFieldName(), vtkm::cont::Field::Association::WholeDataSet, binArray });
|
|
|
|
|
|
|
|
// The output is a "summary" of the input, no need to map fields
|
2018-04-11 00:11:32 +00:00
|
|
|
return output;
|
|
|
|
}
|
|
|
|
|
2022-09-26 16:23:26 +00:00
|
|
|
VTKM_CONT vtkm::cont::PartitionedDataSet HistogramMPI::DoExecutePartitions(
|
|
|
|
const vtkm::cont::PartitionedDataSet& input)
|
|
|
|
{
|
|
|
|
this->PreExecute(input);
|
2022-12-01 20:07:56 +00:00
|
|
|
auto result = this->Filter::DoExecutePartitions(input);
|
2022-09-26 16:23:26 +00:00
|
|
|
this->PostExecute(input, result);
|
|
|
|
return result;
|
|
|
|
}
|
|
|
|
|
2018-04-11 00:11:32 +00:00
|
|
|
//-----------------------------------------------------------------------------
|
2022-09-26 16:23:26 +00:00
|
|
|
inline VTKM_CONT void HistogramMPI::PreExecute(const vtkm::cont::PartitionedDataSet& input)
|
2018-04-11 00:11:32 +00:00
|
|
|
{
|
|
|
|
if (this->Range.IsNonEmpty())
|
|
|
|
{
|
|
|
|
this->ComputedRange = this->Range;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
auto handle = vtkm::cont::FieldRangeGlobalCompute(
|
|
|
|
input, this->GetActiveFieldName(), this->GetActiveFieldAssociation());
|
|
|
|
if (handle.GetNumberOfValues() != 1)
|
|
|
|
{
|
|
|
|
throw vtkm::cont::ErrorFilterExecution("expecting scalar field.");
|
|
|
|
}
|
2020-01-28 19:14:32 +00:00
|
|
|
this->ComputedRange = handle.ReadPortal().Get(0);
|
2018-04-11 00:11:32 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
//-----------------------------------------------------------------------------
|
2019-09-02 14:38:47 +00:00
|
|
|
inline VTKM_CONT void HistogramMPI::PostExecute(const vtkm::cont::PartitionedDataSet&,
|
2022-09-26 16:23:26 +00:00
|
|
|
vtkm::cont::PartitionedDataSet& result)
|
2018-04-11 00:11:32 +00:00
|
|
|
{
|
|
|
|
// iterate and compute HistogramMPI for each local block.
|
2019-09-02 14:38:47 +00:00
|
|
|
detail::DistributedHistogram helper(result.GetNumberOfPartitions());
|
|
|
|
for (vtkm::Id cc = 0; cc < result.GetNumberOfPartitions(); ++cc)
|
2018-04-11 00:11:32 +00:00
|
|
|
{
|
2019-09-02 14:38:47 +00:00
|
|
|
auto& ablock = result.GetPartition(cc);
|
2018-04-11 00:11:32 +00:00
|
|
|
helper.SetLocalHistogram(cc, ablock.GetField(this->GetOutputFieldName()));
|
|
|
|
}
|
|
|
|
|
|
|
|
vtkm::cont::DataSet output;
|
|
|
|
vtkm::cont::Field rfield(this->GetOutputFieldName(),
|
2022-09-26 16:23:26 +00:00
|
|
|
vtkm::cont::Field::Association::WholeDataSet,
|
2018-04-11 00:11:32 +00:00
|
|
|
helper.ReduceAll(this->NumberOfBins));
|
|
|
|
output.AddField(rfield);
|
|
|
|
|
2019-09-02 14:38:47 +00:00
|
|
|
result = vtkm::cont::PartitionedDataSet(output);
|
2018-04-11 00:11:32 +00:00
|
|
|
}
|
|
|
|
} // namespace example
|