2014-02-11 17:34:56 +00:00
|
|
|
//============================================================================
|
|
|
|
// 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.
|
|
|
|
//
|
|
|
|
// Copyright 2014 Sandia Corporation.
|
|
|
|
// Copyright 2014 UT-Battelle, LLC.
|
2015-05-21 12:09:22 +00:00
|
|
|
// Copyright 2014 Los Alamos National Security.
|
2014-02-11 17:34:56 +00:00
|
|
|
//
|
|
|
|
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
|
|
|
|
// the U.S. Government retains certain rights in this software.
|
|
|
|
//
|
|
|
|
// Under the terms of Contract DE-AC52-06NA25396 with Los Alamos National
|
|
|
|
// Laboratory (LANL), the U.S. Government retains certain rights in
|
|
|
|
// this software.
|
|
|
|
//============================================================================
|
2014-03-07 15:19:09 +00:00
|
|
|
#ifndef vtk_m_cont_internal_DeviceAdapterAlgorithmGeneral_h
|
|
|
|
#define vtk_m_cont_internal_DeviceAdapterAlgorithmGeneral_h
|
2014-02-11 17:34:56 +00:00
|
|
|
|
|
|
|
#include <vtkm/cont/ArrayHandle.h>
|
2015-04-28 13:31:50 +00:00
|
|
|
#include <vtkm/cont/ArrayHandleImplicit.h>
|
2015-09-15 04:11:09 +00:00
|
|
|
#include <vtkm/cont/ArrayHandleIndex.h>
|
2015-04-23 17:25:37 +00:00
|
|
|
#include <vtkm/cont/ArrayHandleZip.h>
|
2015-07-16 19:40:22 +00:00
|
|
|
#include <vtkm/cont/internal/FunctorsGeneral.h>
|
2014-06-10 17:35:13 +00:00
|
|
|
|
2014-02-11 17:34:56 +00:00
|
|
|
#include <vtkm/exec/internal/ErrorMessageBuffer.h>
|
|
|
|
|
2016-03-14 13:51:17 +00:00
|
|
|
VTKM_THIRDPARTY_PRE_INCLUDE
|
|
|
|
#if defined(VTKM_MSVC)
|
|
|
|
#define WIN32_LEAN_AND_MEAN
|
|
|
|
#define NOMINMAX
|
|
|
|
#include <Windows.h>
|
|
|
|
#undef WIN32_LEAN_AND_MEAN
|
|
|
|
#undef NOMINMAX
|
|
|
|
#endif
|
|
|
|
VTKM_THIRDPARTY_POST_INCLUDE
|
|
|
|
|
2014-02-11 17:34:56 +00:00
|
|
|
namespace vtkm {
|
|
|
|
namespace cont {
|
|
|
|
namespace internal {
|
|
|
|
|
|
|
|
/// \brief
|
|
|
|
///
|
|
|
|
/// This struct provides algorithms that implement "general" device adapter
|
|
|
|
/// algorithms. If a device adapter provides implementations for Schedule,
|
|
|
|
/// and Synchronize, the rest of the algorithms can be implemented by calling
|
|
|
|
/// these functions.
|
|
|
|
///
|
|
|
|
/// It should be noted that we recommend that you also implement Sort,
|
|
|
|
/// ScanInclusive, and ScanExclusive for improved performance.
|
|
|
|
///
|
|
|
|
/// An easy way to implement the DeviceAdapterAlgorithm specialization is to
|
|
|
|
/// subclass this and override the implementation of methods as necessary.
|
|
|
|
/// As an example, the code would look something like this.
|
|
|
|
///
|
|
|
|
/// \code{.cpp}
|
|
|
|
/// template<>
|
|
|
|
/// struct DeviceAdapterAlgorithm<DeviceAdapterTagFoo>
|
|
|
|
/// : DeviceAdapterAlgorithmGeneral<DeviceAdapterAlgorithm<DeviceAdapterTagFoo>,
|
|
|
|
/// DeviceAdapterTagFoo>
|
|
|
|
/// {
|
|
|
|
/// template<class Functor>
|
|
|
|
/// VTKM_CONT_EXPORT static void Schedule(Functor functor,
|
|
|
|
/// vtkm::Id numInstances)
|
|
|
|
/// {
|
|
|
|
/// ...
|
|
|
|
/// }
|
|
|
|
///
|
|
|
|
/// template<class Functor>
|
|
|
|
/// VTKM_CONT_EXPORT static void Schedule(Functor functor,
|
|
|
|
/// vtkm::Id3 maxRange)
|
|
|
|
/// {
|
|
|
|
/// ...
|
|
|
|
/// }
|
|
|
|
///
|
|
|
|
/// VTKM_CONT_EXPORT static void Synchronize()
|
|
|
|
/// {
|
|
|
|
/// ...
|
|
|
|
/// }
|
|
|
|
/// };
|
|
|
|
/// \endcode
|
|
|
|
///
|
|
|
|
/// You might note that DeviceAdapterAlgorithmGeneral has two template
|
|
|
|
/// parameters that are redundant. Although the first parameter, the class for
|
|
|
|
/// the actual DeviceAdapterAlgorithm class containing Schedule, and
|
|
|
|
/// Synchronize is the same as DeviceAdapterAlgorithm<DeviceAdapterTag>, it is
|
|
|
|
/// made a separate template parameter to avoid a recursive dependence between
|
|
|
|
/// DeviceAdapterAlgorithmGeneral.h and DeviceAdapterAlgorithm.h
|
|
|
|
///
|
|
|
|
template<class DerivedAlgorithm, class DeviceAdapterTag>
|
|
|
|
struct DeviceAdapterAlgorithmGeneral
|
|
|
|
{
|
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Get Execution Value
|
|
|
|
// This method is used internally to get a single element from the execution
|
|
|
|
// array. Might want to expose this and/or allow actual device adapter
|
|
|
|
// implementations to provide one.
|
|
|
|
private:
|
|
|
|
template<typename T, class CIn>
|
|
|
|
VTKM_CONT_EXPORT
|
|
|
|
static T GetExecutionValue(const vtkm::cont::ArrayHandle<T, CIn> &input,
|
|
|
|
vtkm::Id index)
|
|
|
|
{
|
|
|
|
typedef vtkm::cont::ArrayHandle<T,CIn> InputArrayType;
|
2014-06-23 23:33:04 +00:00
|
|
|
typedef vtkm::cont::ArrayHandle<T,vtkm::cont::StorageTagBasic>
|
2014-02-11 17:34:56 +00:00
|
|
|
OutputArrayType;
|
|
|
|
|
|
|
|
OutputArrayType output;
|
|
|
|
|
|
|
|
CopyKernel<
|
|
|
|
typename InputArrayType::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
|
|
|
typename OutputArrayType::template ExecutionTypes<DeviceAdapterTag>::Portal>
|
|
|
|
kernel(input.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
output.PrepareForOutput(1, DeviceAdapterTag()),
|
|
|
|
index);
|
|
|
|
|
|
|
|
DerivedAlgorithm::Schedule(kernel, 1);
|
|
|
|
|
|
|
|
return output.GetPortalConstControl().Get(0);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
public:
|
2015-07-16 19:40:22 +00:00
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Copy
|
2015-06-12 15:46:36 +00:00
|
|
|
template<typename T, typename U, class CIn, class COut>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static void Copy(const vtkm::cont::ArrayHandle<T, CIn> &input,
|
2015-06-12 15:46:36 +00:00
|
|
|
vtkm::cont::ArrayHandle<U, COut> &output)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
|
|
|
vtkm::Id arraySize = input.GetNumberOfValues();
|
|
|
|
|
|
|
|
CopyKernel<
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CIn>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
2015-06-12 15:46:36 +00:00
|
|
|
typename vtkm::cont::ArrayHandle<U,COut>::template ExecutionTypes<DeviceAdapterTag>::Portal>
|
2014-02-11 17:34:56 +00:00
|
|
|
kernel(input.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
output.PrepareForOutput(arraySize, DeviceAdapterTag()));
|
|
|
|
|
|
|
|
DerivedAlgorithm::Schedule(kernel, arraySize);
|
|
|
|
}
|
|
|
|
|
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Lower Bounds
|
|
|
|
template<typename T, class CIn, class CVal, class COut>
|
|
|
|
VTKM_CONT_EXPORT static void LowerBounds(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
|
|
|
const vtkm::cont::ArrayHandle<T,CVal> &values,
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id,COut> &output)
|
|
|
|
{
|
|
|
|
vtkm::Id arraySize = values.GetNumberOfValues();
|
|
|
|
|
|
|
|
LowerBoundsKernel<
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CIn>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CVal>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
|
|
|
typename vtkm::cont::ArrayHandle<vtkm::Id,COut>::template ExecutionTypes<DeviceAdapterTag>::Portal>
|
|
|
|
kernel(input.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
values.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
output.PrepareForOutput(arraySize, DeviceAdapterTag()));
|
|
|
|
|
|
|
|
DerivedAlgorithm::Schedule(kernel, arraySize);
|
|
|
|
}
|
|
|
|
|
2015-06-22 15:49:11 +00:00
|
|
|
template<typename T, class CIn, class CVal, class COut, class BinaryCompare>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static void LowerBounds(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
|
|
|
const vtkm::cont::ArrayHandle<T,CVal> &values,
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id,COut> &output,
|
2015-06-22 15:49:11 +00:00
|
|
|
BinaryCompare binary_compare)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
|
|
|
vtkm::Id arraySize = values.GetNumberOfValues();
|
|
|
|
|
|
|
|
LowerBoundsComparisonKernel<
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CIn>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CVal>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
|
|
|
typename vtkm::cont::ArrayHandle<vtkm::Id,COut>::template ExecutionTypes<DeviceAdapterTag>::Portal,
|
2015-06-22 15:49:11 +00:00
|
|
|
BinaryCompare>
|
2014-02-11 17:34:56 +00:00
|
|
|
kernel(input.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
values.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
output.PrepareForOutput(arraySize, DeviceAdapterTag()),
|
2015-06-22 15:49:11 +00:00
|
|
|
binary_compare);
|
2014-02-11 17:34:56 +00:00
|
|
|
|
|
|
|
DerivedAlgorithm::Schedule(kernel, arraySize);
|
|
|
|
}
|
|
|
|
|
|
|
|
template<class CIn, class COut>
|
|
|
|
VTKM_CONT_EXPORT static void LowerBounds(
|
|
|
|
const vtkm::cont::ArrayHandle<vtkm::Id,CIn> &input,
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id,COut> &values_output)
|
|
|
|
{
|
|
|
|
DeviceAdapterAlgorithmGeneral<
|
|
|
|
DerivedAlgorithm,DeviceAdapterTag>::LowerBounds(input,
|
|
|
|
values_output,
|
|
|
|
values_output);
|
|
|
|
}
|
|
|
|
|
2015-04-28 13:31:50 +00:00
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Reduce
|
|
|
|
template<typename T, class CIn>
|
|
|
|
VTKM_CONT_EXPORT static T Reduce(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input, T initialValue)
|
2015-04-30 13:01:52 +00:00
|
|
|
{
|
|
|
|
return DerivedAlgorithm::Reduce(input, initialValue, vtkm::internal::Add());
|
|
|
|
}
|
|
|
|
|
2015-06-22 15:55:02 +00:00
|
|
|
template<typename T, class CIn, class BinaryFunctor>
|
2015-04-30 13:01:52 +00:00
|
|
|
VTKM_CONT_EXPORT static T Reduce(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
|
|
|
T initialValue,
|
2015-06-22 15:55:02 +00:00
|
|
|
BinaryFunctor binary_functor)
|
2015-04-28 13:31:50 +00:00
|
|
|
{
|
|
|
|
//Crazy Idea:
|
|
|
|
//We create a implicit array handle that wraps the input
|
|
|
|
//array handle. The implicit functor is passed the input array handle, and
|
|
|
|
//the number of elements it needs to sum. This way the implicit handle
|
|
|
|
//acts as the first level reduction. Say for example reducing 16 values
|
|
|
|
//at a time.
|
|
|
|
//
|
|
|
|
//Now that we have an implicit array that is 1/16 the length of full array
|
|
|
|
//we can use scan inclusive to compute the final sum
|
2015-07-16 19:40:22 +00:00
|
|
|
typedef typename vtkm::cont::ArrayHandle<T,CIn>::template ExecutionTypes<DeviceAdapterTag>
|
|
|
|
::PortalConst InputPortalType;
|
|
|
|
|
2015-04-28 13:31:50 +00:00
|
|
|
typedef ReduceKernel<
|
2015-07-16 19:40:22 +00:00
|
|
|
InputPortalType,
|
2015-06-22 15:55:02 +00:00
|
|
|
BinaryFunctor
|
2015-04-28 13:31:50 +00:00
|
|
|
> ReduceKernelType;
|
|
|
|
|
|
|
|
typedef vtkm::cont::ArrayHandleImplicit<
|
|
|
|
T,
|
|
|
|
ReduceKernelType > ReduceHandleType;
|
|
|
|
typedef vtkm::cont::ArrayHandle<
|
|
|
|
T,
|
|
|
|
vtkm::cont::StorageTagBasic> TempArrayType;
|
|
|
|
|
2015-07-16 19:40:22 +00:00
|
|
|
ReduceKernelType kernel(input.PrepareForInput( DeviceAdapterTag() ),
|
|
|
|
binary_functor);
|
|
|
|
|
2015-04-28 13:31:50 +00:00
|
|
|
vtkm::Id length = (input.GetNumberOfValues() / 16);
|
|
|
|
length += (input.GetNumberOfValues() % 16 == 0) ? 0 : 1;
|
|
|
|
ReduceHandleType reduced = vtkm::cont::make_ArrayHandleImplicit<T>(kernel,
|
|
|
|
length);
|
|
|
|
|
2015-04-30 13:01:52 +00:00
|
|
|
TempArrayType inclusiveScanStorage;
|
|
|
|
T scanResult = DerivedAlgorithm::ScanInclusive(reduced,
|
|
|
|
inclusiveScanStorage,
|
2015-06-22 15:55:02 +00:00
|
|
|
binary_functor);
|
|
|
|
return binary_functor(initialValue, scanResult);
|
2015-04-28 13:31:50 +00:00
|
|
|
}
|
|
|
|
|
2015-05-04 19:53:35 +00:00
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Reduce By Key
|
|
|
|
template<typename T, typename U, class KIn, class VIn, class KOut, class VOut,
|
2015-06-22 15:55:02 +00:00
|
|
|
class BinaryFunctor>
|
2015-05-04 19:53:35 +00:00
|
|
|
VTKM_CONT_EXPORT static void ReduceByKey(
|
|
|
|
const vtkm::cont::ArrayHandle<T,KIn> &keys,
|
|
|
|
const vtkm::cont::ArrayHandle<U,VIn> &values,
|
|
|
|
vtkm::cont::ArrayHandle<T,KOut> &keys_output,
|
|
|
|
vtkm::cont::ArrayHandle<U,VOut> &values_output,
|
2015-06-22 15:55:02 +00:00
|
|
|
BinaryFunctor binary_functor)
|
2015-05-04 19:53:35 +00:00
|
|
|
{
|
2016-04-20 21:41:14 +00:00
|
|
|
VTKM_ASSERT(keys.GetNumberOfValues() == values.GetNumberOfValues());
|
2015-05-04 19:53:35 +00:00
|
|
|
const vtkm::Id numberOfKeys = keys.GetNumberOfValues();
|
|
|
|
|
|
|
|
if(numberOfKeys <= 1)
|
|
|
|
{ //we only have a single key/value so that is our output
|
|
|
|
DerivedAlgorithm::Copy(keys, keys_output);
|
|
|
|
DerivedAlgorithm::Copy(values, values_output);
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
//we need to determine based on the keys what is the keystate for
|
|
|
|
//each key. The states are start, middle, end of a series and the special
|
|
|
|
//state start and end of a series
|
2015-06-12 20:38:03 +00:00
|
|
|
vtkm::cont::ArrayHandle< ReduceKeySeriesStates > keystate;
|
2015-05-04 19:53:35 +00:00
|
|
|
|
|
|
|
{
|
|
|
|
typedef typename vtkm::cont::ArrayHandle<T,KIn>::template ExecutionTypes<DeviceAdapterTag>
|
|
|
|
::PortalConst InputPortalType;
|
|
|
|
|
2015-06-12 20:38:03 +00:00
|
|
|
typedef typename vtkm::cont::ArrayHandle< ReduceKeySeriesStates >::template ExecutionTypes<DeviceAdapterTag>
|
2015-05-04 19:53:35 +00:00
|
|
|
::Portal KeyStatePortalType;
|
|
|
|
|
|
|
|
InputPortalType inputPortal = keys.PrepareForInput(DeviceAdapterTag());
|
|
|
|
KeyStatePortalType keyStatePortal = keystate.PrepareForOutput(numberOfKeys,
|
|
|
|
DeviceAdapterTag());
|
2015-07-16 19:40:22 +00:00
|
|
|
ReduceStencilGeneration<InputPortalType, KeyStatePortalType> kernel(inputPortal, keyStatePortal);
|
2015-05-04 19:53:35 +00:00
|
|
|
DerivedAlgorithm::Schedule(kernel, numberOfKeys);
|
|
|
|
}
|
|
|
|
|
|
|
|
//next step is we need to reduce the values for each key. This is done
|
|
|
|
//by running an inclusive scan over the values array using the stencil.
|
|
|
|
//
|
|
|
|
// this inclusive scan will write out two values, the first being
|
|
|
|
// the value summed currently, the second being 0 or 1, with 1 being used
|
|
|
|
// when this is a value of a key we need to write ( END or START_AND_END)
|
|
|
|
{
|
2015-06-05 18:12:07 +00:00
|
|
|
typedef vtkm::cont::ArrayHandle<U,VIn> ValueInHandleType;
|
|
|
|
typedef vtkm::cont::ArrayHandle<U,VOut> ValueOutHandleType;
|
2015-06-12 20:38:03 +00:00
|
|
|
typedef vtkm::cont::ArrayHandle< ReduceKeySeriesStates> StencilHandleType;
|
2015-06-05 18:12:07 +00:00
|
|
|
typedef vtkm::cont::ArrayHandleZip<ValueInHandleType,
|
|
|
|
StencilHandleType> ZipInHandleType;
|
|
|
|
typedef vtkm::cont::ArrayHandleZip<ValueOutHandleType,
|
|
|
|
StencilHandleType> ZipOutHandleType;
|
2015-05-04 19:53:35 +00:00
|
|
|
|
2015-06-05 18:12:07 +00:00
|
|
|
StencilHandleType stencil;
|
|
|
|
ValueOutHandleType reducedValues;
|
2015-05-04 19:53:35 +00:00
|
|
|
|
2015-06-05 18:12:07 +00:00
|
|
|
ZipInHandleType scanInput( values, keystate);
|
|
|
|
ZipOutHandleType scanOutput( reducedValues, stencil);
|
2015-06-12 20:11:28 +00:00
|
|
|
|
|
|
|
DerivedAlgorithm::ScanInclusive(scanInput,
|
|
|
|
scanOutput,
|
2015-06-22 15:55:02 +00:00
|
|
|
ReduceByKeyAdd<BinaryFunctor>(binary_functor) );
|
2015-05-04 19:53:35 +00:00
|
|
|
|
|
|
|
//at this point we are done with keystate, so free the memory
|
|
|
|
keystate.ReleaseResources();
|
|
|
|
|
|
|
|
// all we need know is an efficient way of doing the write back to the
|
|
|
|
// reduced global memory. this is done by using StreamCompact with the
|
|
|
|
// stencil and values we just created with the inclusive scan
|
|
|
|
DerivedAlgorithm::StreamCompact( reducedValues,
|
|
|
|
stencil,
|
|
|
|
values_output,
|
|
|
|
ReduceByKeyUnaryStencilOp());
|
|
|
|
|
|
|
|
} //release all temporary memory
|
|
|
|
|
|
|
|
|
|
|
|
//find all the unique keys
|
|
|
|
DerivedAlgorithm::Copy(keys,keys_output);
|
|
|
|
DerivedAlgorithm::Unique(keys_output);
|
|
|
|
}
|
|
|
|
|
2014-02-11 17:34:56 +00:00
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Scan Exclusive
|
2015-08-27 17:50:32 +00:00
|
|
|
template<typename T, class CIn, class COut, class BinaryFunctor>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static T ScanExclusive(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
2015-08-27 17:50:32 +00:00
|
|
|
vtkm::cont::ArrayHandle<T,COut>& output,
|
|
|
|
BinaryFunctor binaryFunctor,
|
|
|
|
const T& initialValue)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
2014-06-23 23:33:04 +00:00
|
|
|
typedef vtkm::cont::ArrayHandle<T,vtkm::cont::StorageTagBasic>
|
2014-02-11 17:34:56 +00:00
|
|
|
TempArrayType;
|
|
|
|
typedef vtkm::cont::ArrayHandle<T,COut> OutputArrayType;
|
|
|
|
|
2015-08-27 17:50:32 +00:00
|
|
|
typedef typename TempArrayType::template ExecutionTypes<DeviceAdapterTag>
|
|
|
|
::PortalConst SrcPortalType;
|
|
|
|
typedef typename OutputArrayType::template ExecutionTypes<DeviceAdapterTag>
|
|
|
|
::Portal DestPortalType;
|
2014-02-11 17:34:56 +00:00
|
|
|
|
2015-08-27 17:50:32 +00:00
|
|
|
vtkm::Id numValues = input.GetNumberOfValues();
|
|
|
|
if (numValues <= 0)
|
2014-02-11 21:20:30 +00:00
|
|
|
{
|
2015-08-27 17:50:32 +00:00
|
|
|
return initialValue;
|
2014-02-11 21:20:30 +00:00
|
|
|
}
|
2014-02-11 17:34:56 +00:00
|
|
|
|
2015-08-27 17:50:32 +00:00
|
|
|
TempArrayType inclusiveScan;
|
|
|
|
T result = DerivedAlgorithm::ScanInclusive(input, inclusiveScan, binaryFunctor);
|
2014-02-11 17:34:56 +00:00
|
|
|
|
2015-08-27 17:50:32 +00:00
|
|
|
InclusiveToExclusiveKernel<SrcPortalType, DestPortalType, BinaryFunctor>
|
|
|
|
inclusiveToExclusive(inclusiveScan.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
output.PrepareForOutput(numValues, DeviceAdapterTag()),
|
|
|
|
binaryFunctor,
|
|
|
|
initialValue);
|
|
|
|
|
|
|
|
DerivedAlgorithm::Schedule(inclusiveToExclusive, numValues);
|
|
|
|
|
|
|
|
return binaryFunctor(initialValue, result);
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename T, class CIn, class COut>
|
|
|
|
VTKM_CONT_EXPORT static T ScanExclusive(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
|
|
|
vtkm::cont::ArrayHandle<T,COut>& output)
|
|
|
|
{
|
|
|
|
return ScanExclusive(input, output, vtkm::Sum(),
|
|
|
|
vtkm::TypeTraits<T>::ZeroInitialization());
|
2014-02-11 17:34:56 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Scan Inclusive
|
|
|
|
template<typename T, class CIn, class COut>
|
|
|
|
VTKM_CONT_EXPORT static T ScanInclusive(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
|
|
|
vtkm::cont::ArrayHandle<T,COut>& output)
|
2015-04-29 19:48:25 +00:00
|
|
|
{
|
|
|
|
return DerivedAlgorithm::ScanInclusive(input,
|
|
|
|
output,
|
|
|
|
vtkm::internal::Add());
|
|
|
|
}
|
|
|
|
|
2015-06-23 12:25:06 +00:00
|
|
|
template<typename T, class CIn, class COut, class BinaryFunctor>
|
2015-04-29 19:48:25 +00:00
|
|
|
VTKM_CONT_EXPORT static T ScanInclusive(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
|
|
|
vtkm::cont::ArrayHandle<T,COut>& output,
|
2015-06-23 12:25:06 +00:00
|
|
|
BinaryFunctor binary_functor)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
|
|
|
typedef typename
|
|
|
|
vtkm::cont::ArrayHandle<T,COut>
|
|
|
|
::template ExecutionTypes<DeviceAdapterTag>::Portal PortalType;
|
|
|
|
|
2015-06-23 12:25:06 +00:00
|
|
|
typedef ScanKernel<PortalType,BinaryFunctor> ScanKernelType;
|
2015-04-29 19:48:25 +00:00
|
|
|
|
2014-02-11 17:34:56 +00:00
|
|
|
DerivedAlgorithm::Copy(input, output);
|
|
|
|
|
|
|
|
vtkm::Id numValues = output.GetNumberOfValues();
|
|
|
|
if (numValues < 1)
|
2015-05-20 13:29:31 +00:00
|
|
|
{
|
|
|
|
return output.GetPortalConstControl().Get(0);
|
|
|
|
}
|
2014-02-11 17:34:56 +00:00
|
|
|
|
|
|
|
PortalType portal = output.PrepareForInPlace(DeviceAdapterTag());
|
|
|
|
|
|
|
|
vtkm::Id stride;
|
|
|
|
for (stride = 2; stride-1 < numValues; stride *= 2)
|
2014-02-11 21:20:30 +00:00
|
|
|
{
|
2015-06-23 12:25:06 +00:00
|
|
|
ScanKernelType kernel(portal, binary_functor, stride, stride/2 - 1);
|
2014-02-11 17:34:56 +00:00
|
|
|
DerivedAlgorithm::Schedule(kernel, numValues/stride);
|
2014-02-11 21:20:30 +00:00
|
|
|
}
|
2014-02-11 17:34:56 +00:00
|
|
|
|
|
|
|
// Do reverse operation on odd indices. Start at stride we were just at.
|
|
|
|
for (stride /= 2; stride > 1; stride /= 2)
|
2014-02-11 21:20:30 +00:00
|
|
|
{
|
2015-06-23 12:25:06 +00:00
|
|
|
ScanKernelType kernel(portal, binary_functor, stride, stride - 1);
|
2014-02-11 17:34:56 +00:00
|
|
|
DerivedAlgorithm::Schedule(kernel, numValues/stride);
|
2014-02-11 21:20:30 +00:00
|
|
|
}
|
2014-02-11 17:34:56 +00:00
|
|
|
|
|
|
|
return GetExecutionValue(output, numValues-1);
|
|
|
|
}
|
|
|
|
|
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Sort
|
2015-06-23 12:33:02 +00:00
|
|
|
template<typename T, class Storage, class BinaryCompare>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static void Sort(
|
2014-06-23 23:33:04 +00:00
|
|
|
vtkm::cont::ArrayHandle<T,Storage> &values,
|
2015-06-23 12:33:02 +00:00
|
|
|
BinaryCompare binary_compare)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
2014-06-23 23:33:04 +00:00
|
|
|
typedef typename vtkm::cont::ArrayHandle<T,Storage> ArrayType;
|
2014-02-11 17:34:56 +00:00
|
|
|
typedef typename ArrayType::template ExecutionTypes<DeviceAdapterTag>
|
|
|
|
::Portal PortalType;
|
|
|
|
|
|
|
|
vtkm::Id numValues = values.GetNumberOfValues();
|
|
|
|
if (numValues < 2) { return; }
|
|
|
|
|
|
|
|
PortalType portal = values.PrepareForInPlace(DeviceAdapterTag());
|
|
|
|
|
|
|
|
vtkm::Id numThreads = 1;
|
|
|
|
while (numThreads < numValues) { numThreads *= 2; }
|
|
|
|
numThreads /= 2;
|
|
|
|
|
2015-06-23 12:33:02 +00:00
|
|
|
typedef BitonicSortMergeKernel<PortalType,BinaryCompare> MergeKernel;
|
|
|
|
typedef BitonicSortCrossoverKernel<PortalType,BinaryCompare> CrossoverKernel;
|
2014-02-11 17:34:56 +00:00
|
|
|
|
|
|
|
for (vtkm::Id crossoverSize = 1;
|
|
|
|
crossoverSize < numValues;
|
|
|
|
crossoverSize *= 2)
|
2014-02-11 21:20:30 +00:00
|
|
|
{
|
2015-06-23 12:33:02 +00:00
|
|
|
DerivedAlgorithm::Schedule(CrossoverKernel(portal,binary_compare,crossoverSize),
|
2014-02-11 17:34:56 +00:00
|
|
|
numThreads);
|
|
|
|
for (vtkm::Id mergeSize = crossoverSize/2; mergeSize > 0; mergeSize /= 2)
|
2014-02-11 21:20:30 +00:00
|
|
|
{
|
2015-06-23 12:33:02 +00:00
|
|
|
DerivedAlgorithm::Schedule(MergeKernel(portal,binary_compare,mergeSize),
|
2014-02-11 17:34:56 +00:00
|
|
|
numThreads);
|
|
|
|
}
|
2014-02-11 21:20:30 +00:00
|
|
|
}
|
2014-02-11 17:34:56 +00:00
|
|
|
}
|
|
|
|
|
2014-06-23 23:33:04 +00:00
|
|
|
template<typename T, class Storage>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static void Sort(
|
2014-06-23 23:33:04 +00:00
|
|
|
vtkm::cont::ArrayHandle<T,Storage> &values)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
|
|
|
DerivedAlgorithm::Sort(values, DefaultCompareFunctor());
|
2015-04-23 17:25:37 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Sort by Key
|
|
|
|
public:
|
|
|
|
|
|
|
|
template<typename T, typename U, class StorageT, class StorageU>
|
|
|
|
VTKM_CONT_EXPORT static void SortByKey(
|
|
|
|
vtkm::cont::ArrayHandle<T,StorageT> &keys,
|
|
|
|
vtkm::cont::ArrayHandle<U,StorageU> &values)
|
|
|
|
{
|
|
|
|
//combine the keys and values into a ZipArrayHandle
|
|
|
|
//we than need to specify a custom compare function wrapper
|
|
|
|
//that only checks for key side of the pair, using a custom compare functor.
|
|
|
|
typedef vtkm::cont::ArrayHandle<T,StorageT> KeyType;
|
|
|
|
typedef vtkm::cont::ArrayHandle<U,StorageU> ValueType;
|
|
|
|
typedef vtkm::cont::ArrayHandleZip<KeyType,ValueType> ZipHandleType;
|
|
|
|
|
|
|
|
ZipHandleType zipHandle =
|
|
|
|
vtkm::cont::make_ArrayHandleZip(keys,values);
|
2015-07-16 19:40:22 +00:00
|
|
|
DerivedAlgorithm::Sort(zipHandle,internal::KeyCompare<T,U>());
|
2015-04-23 17:25:37 +00:00
|
|
|
}
|
|
|
|
|
2015-06-23 12:33:02 +00:00
|
|
|
template<typename T, typename U, class StorageT, class StorageU, class BinaryCompare>
|
2015-04-23 17:25:37 +00:00
|
|
|
VTKM_CONT_EXPORT static void SortByKey(
|
|
|
|
vtkm::cont::ArrayHandle<T,StorageT> &keys,
|
|
|
|
vtkm::cont::ArrayHandle<U,StorageU> &values,
|
2015-06-23 12:33:02 +00:00
|
|
|
BinaryCompare binary_compare)
|
2015-04-23 17:25:37 +00:00
|
|
|
{
|
|
|
|
//combine the keys and values into a ZipArrayHandle
|
|
|
|
//we than need to specify a custom compare function wrapper
|
|
|
|
//that only checks for key side of the pair, using the custom compare
|
|
|
|
//functor that the user passed in
|
|
|
|
typedef vtkm::cont::ArrayHandle<T,StorageT> KeyType;
|
|
|
|
typedef vtkm::cont::ArrayHandle<U,StorageU> ValueType;
|
|
|
|
typedef vtkm::cont::ArrayHandleZip<KeyType,ValueType> ZipHandleType;
|
|
|
|
|
|
|
|
ZipHandleType zipHandle =
|
|
|
|
vtkm::cont::make_ArrayHandleZip(keys,values);
|
2015-07-16 19:40:22 +00:00
|
|
|
DerivedAlgorithm::Sort(zipHandle,internal::KeyCompare<T,U,BinaryCompare>(binary_compare));
|
2014-02-11 17:34:56 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Stream Compact
|
2015-05-20 15:13:42 +00:00
|
|
|
template<typename T, typename U, class CIn, class CStencil,
|
2015-06-23 12:56:18 +00:00
|
|
|
class COut, class UnaryPredicate>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static void StreamCompact(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn>& input,
|
|
|
|
const vtkm::cont::ArrayHandle<U,CStencil>& stencil,
|
2015-05-20 15:13:42 +00:00
|
|
|
vtkm::cont::ArrayHandle<T,COut>& output,
|
2015-06-23 12:56:18 +00:00
|
|
|
UnaryPredicate unary_predicate)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
2016-04-20 21:41:14 +00:00
|
|
|
VTKM_ASSERT(input.GetNumberOfValues() == stencil.GetNumberOfValues());
|
2014-02-11 17:34:56 +00:00
|
|
|
vtkm::Id arrayLength = stencil.GetNumberOfValues();
|
|
|
|
|
|
|
|
typedef vtkm::cont::ArrayHandle<
|
2014-06-23 23:33:04 +00:00
|
|
|
vtkm::Id, vtkm::cont::StorageTagBasic> IndexArrayType;
|
2014-02-11 17:34:56 +00:00
|
|
|
IndexArrayType indices;
|
|
|
|
|
|
|
|
typedef typename vtkm::cont::ArrayHandle<U,CStencil>
|
|
|
|
::template ExecutionTypes<DeviceAdapterTag>::PortalConst
|
|
|
|
StencilPortalType;
|
|
|
|
StencilPortalType stencilPortal =
|
|
|
|
stencil.PrepareForInput(DeviceAdapterTag());
|
|
|
|
|
|
|
|
typedef typename IndexArrayType
|
|
|
|
::template ExecutionTypes<DeviceAdapterTag>::Portal IndexPortalType;
|
|
|
|
IndexPortalType indexPortal =
|
|
|
|
indices.PrepareForOutput(arrayLength, DeviceAdapterTag());
|
|
|
|
|
2015-05-20 15:13:42 +00:00
|
|
|
StencilToIndexFlagKernel< StencilPortalType,
|
|
|
|
IndexPortalType,
|
2015-06-23 12:56:18 +00:00
|
|
|
UnaryPredicate> indexKernel(stencilPortal,
|
2015-05-20 15:13:42 +00:00
|
|
|
indexPortal,
|
2015-06-23 12:56:18 +00:00
|
|
|
unary_predicate);
|
2014-02-11 17:34:56 +00:00
|
|
|
|
|
|
|
DerivedAlgorithm::Schedule(indexKernel, arrayLength);
|
|
|
|
|
|
|
|
vtkm::Id outArrayLength = DerivedAlgorithm::ScanExclusive(indices, indices);
|
|
|
|
|
|
|
|
typedef typename vtkm::cont::ArrayHandle<T,CIn>
|
|
|
|
::template ExecutionTypes<DeviceAdapterTag>::PortalConst
|
|
|
|
InputPortalType;
|
|
|
|
InputPortalType inputPortal = input.PrepareForInput(DeviceAdapterTag());
|
|
|
|
|
|
|
|
typedef typename vtkm::cont::ArrayHandle<T,COut>
|
|
|
|
::template ExecutionTypes<DeviceAdapterTag>::Portal OutputPortalType;
|
|
|
|
OutputPortalType outputPortal =
|
|
|
|
output.PrepareForOutput(outArrayLength, DeviceAdapterTag());
|
|
|
|
|
|
|
|
CopyIfKernel<
|
|
|
|
InputPortalType,
|
|
|
|
StencilPortalType,
|
|
|
|
IndexPortalType,
|
2015-05-20 15:13:42 +00:00
|
|
|
OutputPortalType,
|
2015-06-23 12:56:18 +00:00
|
|
|
UnaryPredicate> copyKernel(inputPortal,
|
2014-02-11 17:34:56 +00:00
|
|
|
stencilPortal,
|
|
|
|
indexPortal,
|
2015-05-20 15:13:42 +00:00
|
|
|
outputPortal,
|
2015-06-23 12:56:18 +00:00
|
|
|
unary_predicate);
|
2014-02-11 17:34:56 +00:00
|
|
|
DerivedAlgorithm::Schedule(copyKernel, arrayLength);
|
|
|
|
}
|
|
|
|
|
2015-05-20 15:13:42 +00:00
|
|
|
template<typename T, typename U, class CIn, class CStencil, class COut>
|
|
|
|
VTKM_CONT_EXPORT static void StreamCompact(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn>& input,
|
|
|
|
const vtkm::cont::ArrayHandle<U,CStencil>& stencil,
|
|
|
|
vtkm::cont::ArrayHandle<T,COut>& output)
|
|
|
|
{
|
2015-07-20 20:29:43 +00:00
|
|
|
::vtkm::NotZeroInitialized unary_predicate;
|
2015-06-23 12:56:18 +00:00
|
|
|
DerivedAlgorithm::StreamCompact(input, stencil, output, unary_predicate);
|
2015-05-20 15:13:42 +00:00
|
|
|
}
|
|
|
|
|
2014-02-11 17:34:56 +00:00
|
|
|
template<typename T, class CStencil, class COut>
|
|
|
|
VTKM_CONT_EXPORT static void StreamCompact(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CStencil> &stencil,
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id,COut> &output)
|
|
|
|
{
|
2015-09-15 04:11:09 +00:00
|
|
|
vtkm::cont::ArrayHandleIndex input(stencil.GetNumberOfValues());
|
2014-02-11 17:34:56 +00:00
|
|
|
DerivedAlgorithm::StreamCompact(input, stencil, output);
|
|
|
|
}
|
|
|
|
|
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Unique
|
2014-06-23 23:33:04 +00:00
|
|
|
template<typename T, class Storage>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static void Unique(
|
2014-06-23 23:33:04 +00:00
|
|
|
vtkm::cont::ArrayHandle<T,Storage> &values)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
2016-03-18 20:39:21 +00:00
|
|
|
Unique(values, std::equal_to<T>());
|
2014-02-11 17:34:56 +00:00
|
|
|
}
|
|
|
|
|
2015-06-23 12:58:53 +00:00
|
|
|
template<typename T, class Storage, class BinaryCompare>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static void Unique(
|
2014-06-23 23:33:04 +00:00
|
|
|
vtkm::cont::ArrayHandle<T,Storage> &values,
|
2015-06-23 12:58:53 +00:00
|
|
|
BinaryCompare binary_compare)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
2014-06-23 23:33:04 +00:00
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id, vtkm::cont::StorageTagBasic>
|
2014-02-11 17:34:56 +00:00
|
|
|
stencilArray;
|
|
|
|
vtkm::Id inputSize = values.GetNumberOfValues();
|
|
|
|
|
2016-03-18 20:39:21 +00:00
|
|
|
typedef internal::WrappedBinaryOperator<bool,BinaryCompare> WrappedBOpType;
|
|
|
|
WrappedBOpType wrappedCompare(binary_compare);
|
|
|
|
|
2014-02-11 17:34:56 +00:00
|
|
|
ClassifyUniqueComparisonKernel<
|
2014-06-23 23:33:04 +00:00
|
|
|
typename vtkm::cont::ArrayHandle<T,Storage>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
|
|
|
typename vtkm::cont::ArrayHandle<vtkm::Id,vtkm::cont::StorageTagBasic>::template ExecutionTypes<DeviceAdapterTag>::Portal,
|
2016-03-18 20:39:21 +00:00
|
|
|
WrappedBOpType>
|
2014-02-11 17:34:56 +00:00
|
|
|
classifyKernel(values.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
stencilArray.PrepareForOutput(inputSize, DeviceAdapterTag()),
|
2016-03-18 20:39:21 +00:00
|
|
|
wrappedCompare);
|
|
|
|
|
2014-02-11 17:34:56 +00:00
|
|
|
DerivedAlgorithm::Schedule(classifyKernel, inputSize);
|
|
|
|
|
2014-06-23 23:33:04 +00:00
|
|
|
vtkm::cont::ArrayHandle<T, vtkm::cont::StorageTagBasic>
|
2014-02-11 17:34:56 +00:00
|
|
|
outputArray;
|
|
|
|
|
|
|
|
DerivedAlgorithm::StreamCompact(values, stencilArray, outputArray);
|
|
|
|
|
|
|
|
DerivedAlgorithm::Copy(outputArray, values);
|
|
|
|
}
|
|
|
|
|
|
|
|
//--------------------------------------------------------------------------
|
|
|
|
// Upper bounds
|
|
|
|
template<typename T, class CIn, class CVal, class COut>
|
|
|
|
VTKM_CONT_EXPORT static void UpperBounds(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
|
|
|
const vtkm::cont::ArrayHandle<T,CVal> &values,
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id,COut> &output)
|
|
|
|
{
|
|
|
|
vtkm::Id arraySize = values.GetNumberOfValues();
|
|
|
|
|
|
|
|
UpperBoundsKernel<
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CIn>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CVal>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
2015-07-08 20:01:30 +00:00
|
|
|
typename vtkm::cont::ArrayHandle<vtkm::Id,COut>::template ExecutionTypes<DeviceAdapterTag>::Portal>
|
2014-02-11 17:34:56 +00:00
|
|
|
kernel(input.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
values.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
output.PrepareForOutput(arraySize, DeviceAdapterTag()));
|
|
|
|
|
|
|
|
DerivedAlgorithm::Schedule(kernel, arraySize);
|
|
|
|
}
|
|
|
|
|
2015-06-23 13:19:09 +00:00
|
|
|
template<typename T, class CIn, class CVal, class COut, class BinaryCompare>
|
2014-02-11 17:34:56 +00:00
|
|
|
VTKM_CONT_EXPORT static void UpperBounds(
|
|
|
|
const vtkm::cont::ArrayHandle<T,CIn> &input,
|
|
|
|
const vtkm::cont::ArrayHandle<T,CVal> &values,
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id,COut> &output,
|
2015-06-23 13:19:09 +00:00
|
|
|
BinaryCompare binary_compare)
|
2014-02-11 17:34:56 +00:00
|
|
|
{
|
|
|
|
vtkm::Id arraySize = values.GetNumberOfValues();
|
|
|
|
|
|
|
|
UpperBoundsKernelComparisonKernel<
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CIn>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
|
|
|
typename vtkm::cont::ArrayHandle<T,CVal>::template ExecutionTypes<DeviceAdapterTag>::PortalConst,
|
2015-07-08 20:01:30 +00:00
|
|
|
typename vtkm::cont::ArrayHandle<vtkm::Id,COut>::template ExecutionTypes<DeviceAdapterTag>::Portal,
|
2015-06-23 13:19:09 +00:00
|
|
|
BinaryCompare>
|
2014-02-11 17:34:56 +00:00
|
|
|
kernel(input.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
values.PrepareForInput(DeviceAdapterTag()),
|
|
|
|
output.PrepareForOutput(arraySize, DeviceAdapterTag()),
|
2015-06-23 13:19:09 +00:00
|
|
|
binary_compare);
|
2014-02-11 17:34:56 +00:00
|
|
|
|
|
|
|
DerivedAlgorithm::Schedule(kernel, arraySize);
|
|
|
|
}
|
|
|
|
|
|
|
|
template<class CIn, class COut>
|
|
|
|
VTKM_CONT_EXPORT static void UpperBounds(
|
|
|
|
const vtkm::cont::ArrayHandle<vtkm::Id,CIn> &input,
|
|
|
|
vtkm::cont::ArrayHandle<vtkm::Id,COut> &values_output)
|
|
|
|
{
|
|
|
|
DeviceAdapterAlgorithmGeneral<DerivedAlgorithm,
|
|
|
|
DeviceAdapterTag>::UpperBounds(input, values_output, values_output);
|
|
|
|
}
|
|
|
|
|
|
|
|
};
|
|
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
}
|
|
|
|
} // namespace vtkm::cont::internal
|
|
|
|
|
2016-02-10 15:51:31 +00:00
|
|
|
namespace vtkm {
|
|
|
|
namespace cont {
|
|
|
|
/// \brief Class providing a device-specific atomic interface.
|
|
|
|
///
|
|
|
|
/// The class provide the actual implementation used by vtkm::exec::AtomicArray.
|
|
|
|
/// A serial default implementation is provided. But each device will have a different
|
|
|
|
/// implementation.
|
|
|
|
///
|
|
|
|
/// Serial requires no form of atomicity
|
|
|
|
///
|
|
|
|
template<typename T, typename DeviceTag>
|
|
|
|
class DeviceAdapterAtomicArrayImplementation
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
VTKM_CONT_EXPORT
|
2016-03-14 13:51:17 +00:00
|
|
|
DeviceAdapterAtomicArrayImplementation(
|
|
|
|
vtkm::cont::ArrayHandle<T, vtkm::cont::StorageTagBasic> handle):
|
|
|
|
Iterators( IteratorsType( handle.PrepareForInPlace(DeviceTag()) ) )
|
2016-02-10 15:51:31 +00:00
|
|
|
{
|
|
|
|
}
|
|
|
|
|
|
|
|
VTKM_EXEC_EXPORT
|
|
|
|
T Add(vtkm::Id index, const T& value) const
|
|
|
|
{
|
2016-03-14 13:51:17 +00:00
|
|
|
T* lockedValue;
|
2016-06-13 20:34:27 +00:00
|
|
|
#if defined(_ITERATOR_DEBUG_LEVEL) && _ITERATOR_DEBUG_LEVEL > 0
|
2016-03-14 13:51:17 +00:00
|
|
|
typedef typename vtkm::cont::ArrayPortalToIterators<PortalType>::IteratorType IteratorType;
|
2016-05-03 22:13:09 +00:00
|
|
|
typename IteratorType::pointer temp =
|
|
|
|
&(*(Iterators.GetBegin() + static_cast<std::ptrdiff_t>(index)));
|
2016-03-14 13:51:17 +00:00
|
|
|
lockedValue = temp;
|
|
|
|
return vtkmAtomicAdd(lockedValue, value);
|
|
|
|
#else
|
|
|
|
lockedValue = (Iterators.GetBegin()+index);
|
|
|
|
return vtkmAtomicAdd(lockedValue, value);
|
|
|
|
#endif
|
2016-02-10 15:51:31 +00:00
|
|
|
}
|
|
|
|
|
2016-03-08 17:41:02 +00:00
|
|
|
VTKM_EXEC_EXPORT
|
|
|
|
T CompareAndSwap(vtkm::Id index, const T& newValue, const T& oldValue) const
|
|
|
|
{
|
2016-03-14 13:51:17 +00:00
|
|
|
T* lockedValue;
|
2016-06-13 20:34:27 +00:00
|
|
|
#if defined(_ITERATOR_DEBUG_LEVEL) && _ITERATOR_DEBUG_LEVEL > 0
|
2016-03-14 13:51:17 +00:00
|
|
|
typedef typename vtkm::cont::ArrayPortalToIterators<PortalType>::IteratorType IteratorType;
|
2016-05-03 22:13:09 +00:00
|
|
|
typename IteratorType::pointer temp =
|
|
|
|
&(*(Iterators.GetBegin()+static_cast<std::ptrdiff_t>(index)));
|
2016-03-14 13:51:17 +00:00
|
|
|
lockedValue = temp;
|
|
|
|
return vtkmCompareAndSwap(lockedValue, newValue, oldValue);
|
|
|
|
#else
|
|
|
|
lockedValue = (Iterators.GetBegin()+index);
|
|
|
|
return vtkmCompareAndSwap(lockedValue, newValue, oldValue);
|
|
|
|
#endif
|
2016-03-08 17:41:02 +00:00
|
|
|
}
|
|
|
|
|
2016-02-10 15:51:31 +00:00
|
|
|
private:
|
2016-03-14 13:51:17 +00:00
|
|
|
typedef typename vtkm::cont::ArrayHandle<T,vtkm::cont::StorageTagBasic>
|
|
|
|
::template ExecutionTypes<DeviceTag>::Portal PortalType;
|
|
|
|
typedef vtkm::cont::ArrayPortalToIterators<PortalType> IteratorsType;
|
|
|
|
IteratorsType Iterators;
|
|
|
|
|
|
|
|
#if defined(VTKM_MSVC) //MSVC atomics
|
|
|
|
VTKM_EXEC_EXPORT
|
|
|
|
vtkm::Int32 vtkmAtomicAdd(vtkm::Int32 *address, const vtkm::Int32 &value) const
|
|
|
|
{
|
|
|
|
return InterlockedExchangeAdd(reinterpret_cast<volatile long *>(address),value);
|
|
|
|
}
|
|
|
|
|
|
|
|
VTKM_EXEC_EXPORT
|
|
|
|
vtkm::Int64 vtkmAtomicAdd(vtkm::Int64 *address, const vtkm::Int64 &value) const
|
|
|
|
{
|
|
|
|
return InterlockedExchangeAdd64(reinterpret_cast<volatile long long *>(address),value);
|
|
|
|
}
|
|
|
|
|
|
|
|
VTKM_EXEC_EXPORT
|
|
|
|
vtkm::Int32 vtkmCompareAndSwap(vtkm::Int32 *address, const vtkm::Int32 &newValue, const vtkm::Int32 &oldValue) const
|
|
|
|
{
|
|
|
|
return InterlockedCompareExchange(reinterpret_cast<volatile long *>(address),newValue,oldValue);
|
|
|
|
}
|
|
|
|
|
|
|
|
VTKM_EXEC_EXPORT
|
|
|
|
vtkm::Int64 vtkmCompareAndSwap(vtkm::Int64 *address,const vtkm::Int64 &newValue, const vtkm::Int64 &oldValue) const
|
|
|
|
{
|
|
|
|
return InterlockedCompareExchange64(reinterpret_cast<volatile long long *>(address),newValue, oldValue);
|
|
|
|
}
|
|
|
|
|
|
|
|
#else //gcc built-in atomics
|
2016-03-08 18:39:23 +00:00
|
|
|
|
|
|
|
VTKM_EXEC_EXPORT
|
2016-03-14 13:51:17 +00:00
|
|
|
vtkm::Int32 vtkmAtomicAdd(vtkm::Int32 *address, const vtkm::Int32 &value) const
|
2016-03-08 18:39:23 +00:00
|
|
|
{
|
2016-03-14 13:51:17 +00:00
|
|
|
return __sync_fetch_and_add(address,value);
|
2016-03-08 18:39:23 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
VTKM_EXEC_EXPORT
|
2016-03-14 13:51:17 +00:00
|
|
|
vtkm::Int64 vtkmAtomicAdd(vtkm::Int64 *address, const vtkm::Int64 &value) const
|
2016-03-08 18:39:23 +00:00
|
|
|
{
|
2016-03-14 13:51:17 +00:00
|
|
|
return __sync_fetch_and_add(address,value);
|
2016-03-08 18:39:23 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
VTKM_EXEC_EXPORT
|
2016-03-14 13:51:17 +00:00
|
|
|
vtkm::Int32 vtkmCompareAndSwap(vtkm::Int32 *address, const vtkm::Int32 &newValue, const vtkm::Int32 &oldValue) const
|
2016-03-08 18:39:23 +00:00
|
|
|
{
|
2016-03-14 13:51:17 +00:00
|
|
|
return __sync_val_compare_and_swap(address,oldValue, newValue);
|
2016-03-08 18:39:23 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
VTKM_EXEC_EXPORT
|
2016-03-14 13:51:17 +00:00
|
|
|
vtkm::Int64 vtkmCompareAndSwap(vtkm::Int64 *address,const vtkm::Int64 &newValue, const vtkm::Int64 &oldValue) const
|
2016-03-08 18:39:23 +00:00
|
|
|
{
|
2016-03-14 13:51:17 +00:00
|
|
|
return __sync_val_compare_and_swap(address,oldValue,newValue);
|
2016-03-08 18:39:23 +00:00
|
|
|
}
|
|
|
|
|
2016-03-14 13:51:17 +00:00
|
|
|
#endif
|
2016-02-10 15:51:31 +00:00
|
|
|
};
|
|
|
|
|
2016-02-10 16:21:38 +00:00
|
|
|
}
|
|
|
|
} // namespace vtkm::cont
|
|
|
|
|
2014-06-11 16:43:36 +00:00
|
|
|
#endif //vtk_m_cont_internal_DeviceAdapterAlgorithmGeneral_h
|