FunctorsOpenMP.h 20.8 KB
Newer Older
1
2
3
4
//============================================================================
//  Copyright (c) Kitware, Inc.
//  All rights reserved.
//  See LICENSE.txt for details.
5
//
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
//  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.
//============================================================================
#ifndef vtk_m_cont_openmp_internal_FunctorsOpenMP_h
#define vtk_m_cont_openmp_internal_FunctorsOpenMP_h

#include <vtkm/cont/openmp/internal/DeviceAdapterTagOpenMP.h>

#include <vtkm/cont/internal/FunctorsGeneral.h>

#include <vtkm/BinaryOperators.h>
#include <vtkm/BinaryPredicates.h>
#include <vtkm/Pair.h>
#include <vtkm/Types.h>
#include <vtkm/cont/ArrayHandle.h>
#include <vtkm/cont/ErrorExecution.h>

#include <omp.h>

#include <algorithm>
#include <type_traits>
#include <vector>

// Wrap all '#pragma omp ...' calls in this macro so we can disable them in
// non-omp builds and avoid a multitude of 'ignoring pragma..." warnings.
#ifdef _OPENMP
33
34
#define VTKM_OPENMP_DIRECTIVE_IMPL(fullDir) _Pragma(#fullDir)
#define VTKM_OPENMP_DIRECTIVE(dir) VTKM_OPENMP_DIRECTIVE_IMPL(omp dir)
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
#else // _OPENMP
#define VTKM_OPENMP_DIRECTIVE(directive)
#endif // _OPENMP

// When defined, supported type / operator combinations will use the OpenMP
// reduction(...) clause. Otherwise, all reductions use the general
// implementation with a manual reduction once the threads complete.
// I don't know how, but the benchmarks currently perform better without the
// specializations.
//#define VTKM_OPENMP_USE_NATIVE_REDUCTION

namespace vtkm
{
namespace cont
{
namespace openmp
{

constexpr static vtkm::Id CACHE_LINE_SIZE = 64;
constexpr static vtkm::Id PAGE_SIZE = 4096;

// Returns ceil(num/den) for integral types
template <typename T>
static constexpr T CeilDivide(const T& numerator, const T& denominator)
{
  return (numerator + denominator - 1) / denominator;
}

// Computes the number of values per chunk. Note that numChunks + chunkSize may
// exceed numVals, so be sure to check upper limits.
static void ComputeChunkSize(const vtkm::Id numVals,
                             const vtkm::Id numThreads,
                             const vtkm::Id chunksPerThread,
                             const vtkm::Id bytesPerValue,
                             vtkm::Id& numChunks,
                             vtkm::Id& valuesPerChunk)
{
  // try to evenly distribute pages across chunks:
  const vtkm::Id bytesIn = numVals * bytesPerValue;
  const vtkm::Id pagesIn = CeilDivide(bytesIn, PAGE_SIZE);
  // If we don't have enough pages to honor chunksPerThread, ignore it:
  numChunks = (pagesIn > numThreads * chunksPerThread) ? numThreads * chunksPerThread : numThreads;
  const vtkm::Id pagesPerChunk = CeilDivide(pagesIn, numChunks);
  valuesPerChunk = CeilDivide(pagesPerChunk * PAGE_SIZE, bytesPerValue);
}

template <typename T, typename U>
82
static void DoCopy(T src, U dst, vtkm::Id numVals, std::true_type)
83
84
85
86
87
88
89
90
91
{
  if (numVals)
  {
    std::copy(src, src + numVals, dst);
  }
}

// Don't use std::copy when type conversion is required because MSVC.
template <typename InIterT, typename OutIterT>
92
static void DoCopy(InIterT inIter, OutIterT outIter, vtkm::Id numVals, std::false_type)
93
94
95
96
97
98
99
100
101
102
{
  using ValueType = typename std::iterator_traits<OutIterT>::value_type;

  for (vtkm::Id i = 0; i < numVals; ++i)
  {
    *(outIter++) = static_cast<ValueType>(*(inIter++));
  }
}

template <typename InIterT, typename OutIterT>
103
static void DoCopy(InIterT inIter, OutIterT outIter, vtkm::Id numVals)
104
105
106
107
108
109
110
111
112
{
  using InValueType = typename std::iterator_traits<InIterT>::value_type;
  using OutValueType = typename std::iterator_traits<OutIterT>::value_type;

  DoCopy(inIter, outIter, numVals, std::is_same<InValueType, OutValueType>());
}


template <typename InPortalT, typename OutPortalT>
113
114
115
116
117
static void CopyHelper(InPortalT inPortal,
                       OutPortalT outPortal,
                       vtkm::Id inStart,
                       vtkm::Id outStart,
                       vtkm::Id numVals)
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
{
  using InValueT = typename InPortalT::ValueType;
  using OutValueT = typename OutPortalT::ValueType;
  constexpr auto isSame = std::is_same<InValueT, OutValueT>();

  auto inIter = vtkm::cont::ArrayPortalToIteratorBegin(inPortal) + inStart;
  auto outIter = vtkm::cont::ArrayPortalToIteratorBegin(outPortal) + outStart;
  vtkm::Id valuesPerChunk;

  VTKM_OPENMP_DIRECTIVE(parallel default(none) shared(inIter, outIter, valuesPerChunk, numVals))
  {

    VTKM_OPENMP_DIRECTIVE(single)
    {
      // Evenly distribute full pages to all threads. We manually chunk the
      // data here so that we can exploit std::copy's memmove optimizations.
      vtkm::Id numChunks;
      ComputeChunkSize(
        numVals, omp_get_num_threads(), 8, sizeof(InValueT), numChunks, valuesPerChunk);
    }

VTKM_OPENMP_DIRECTIVE(for schedule(static))
for (vtkm::Id i = 0; i < numVals; i += valuesPerChunk)
{
  vtkm::Id chunkSize = std::min(numVals - i, valuesPerChunk);
  DoCopy(inIter + i, outIter + i, chunkSize, isSame);
}
  }
}

struct CopyIfHelper
{
  vtkm::Id NumValues;
  vtkm::Id NumThreads;
  vtkm::Id ValueSize;

  vtkm::Id NumChunks;
  vtkm::Id ChunkSize;
  std::vector<vtkm::Id> EndIds;

  CopyIfHelper() = default;

  void Initialize(vtkm::Id numValues, vtkm::Id valueSize)
  {
    this->NumValues = numValues;
163
    this->NumThreads = static_cast<vtkm::Id>(omp_get_num_threads());
164
165
166
167
168
169
170
    this->ValueSize = valueSize;

    // Evenly distribute pages across the threads. We manually chunk the
    // data here so that we can exploit std::copy's memmove optimizations.
    ComputeChunkSize(
      this->NumValues, this->NumThreads, 8, valueSize, this->NumChunks, this->ChunkSize);

171
    this->EndIds.resize(static_cast<std::size_t>(this->NumChunks));
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
  }

  template <typename InIterT, typename StencilIterT, typename OutIterT, typename PredicateT>
  void CopyIf(InIterT inIter,
              StencilIterT stencilIter,
              OutIterT outIter,
              PredicateT pred,
              vtkm::Id chunk)
  {
    vtkm::Id startPos = std::min(chunk * this->ChunkSize, this->NumValues);
    vtkm::Id endPos = std::min((chunk + 1) * this->ChunkSize, this->NumValues);

    vtkm::Id outPos = startPos;
    for (vtkm::Id inPos = startPos; inPos < endPos; ++inPos)
    {
      if (pred(stencilIter[inPos]))
      {
        outIter[outPos++] = inIter[inPos];
      }
    }

193
    this->EndIds[static_cast<std::size_t>(chunk)] = outPos;
194
195
196
197
198
199
200
201
202
  }

  template <typename OutIterT>
  vtkm::Id Reduce(OutIterT data)
  {
    vtkm::Id endPos = this->EndIds.front();
    for (vtkm::Id i = 1; i < this->NumChunks; ++i)
    {
      vtkm::Id chunkStart = std::min(i * this->ChunkSize, this->NumValues);
203
      vtkm::Id chunkEnd = this->EndIds[static_cast<std::size_t>(i)];
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
      vtkm::Id numValuesToCopy = chunkEnd - chunkStart;
      if (numValuesToCopy > 0 && chunkStart != endPos)
      {
        std::copy(data + chunkStart, data + chunkEnd, data + endPos);
      }
      endPos += numValuesToCopy;
    }
    return endPos;
  }
};

#ifdef VTKM_OPENMP_USE_NATIVE_REDUCTION
// OpenMP only declares reduction operations for primitive types. This utility
// detects if a type T is supported.
template <typename T>
struct OpenMPReductionSupported : std::false_type
{
};
template <>
struct OpenMPReductionSupported<Int8> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<UInt8> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<Int16> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<UInt16> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<Int32> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<UInt32> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<Int64> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<UInt64> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<Float32> : std::true_type
{
};
template <>
struct OpenMPReductionSupported<Float64> : std::true_type
{
};
#else
template <typename T>
using OpenMPReductionSupported = std::false_type;
#endif // VTKM_OPENMP_USE_NATIVE_REDUCTION

struct ReduceHelper
{
  // Generic implementation:
  template <typename PortalT, typename ReturnType, typename Functor>
  static ReturnType Execute(PortalT portal, ReturnType init, Functor functorIn, std::false_type)
  {
    internal::WrappedBinaryOperator<ReturnType, Functor> f(functorIn);

    const vtkm::Id numVals = portal.GetNumberOfValues();
    auto data = vtkm::cont::ArrayPortalToIteratorBegin(portal);

    bool doParallel = false;
279
280
    int numThreads = 0;
    std::unique_ptr<ReturnType[]> threadData;
281
282

    VTKM_OPENMP_DIRECTIVE(parallel default(none) firstprivate(f)
283
                            shared(data, doParallel, numThreads, threadData))
284
285
286
287
288
289
    {

      int tid = omp_get_thread_num();

      VTKM_OPENMP_DIRECTIVE(single)
      {
290
        numThreads = omp_get_num_threads();
291
292
293
        if (numVals >= numThreads * 2)
        {
          doParallel = true;
294
          threadData.reset(new ReturnType[numThreads]);
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
        }
      }

      if (doParallel)
      {
        // Use the first (numThreads*2) values for initializing:
        ReturnType accum;
        accum = f(data[2 * tid], data[2 * tid + 1]);

        // Assign each thread chunks of the remaining values for local reduction
        VTKM_OPENMP_DIRECTIVE(for schedule(static))
        for (vtkm::Id i = numThreads * 2; i < numVals; i++)
        {
          accum = f(accum, data[i]);
        }

311
        threadData[static_cast<std::size_t>(tid)] = accum;
312
313
314
315
316
317
      }
    } // end parallel

    if (doParallel)
    {
      // do the final reduction serially:
318
      for (size_t i = 0; i < static_cast<size_t>(numThreads); ++i)
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
      {
        init = f(init, threadData[i]);
      }
    }
    else
    {
      // Not enough threads. Do the entire reduction in serial:
      for (vtkm::Id i = 0; i < numVals; ++i)
      {
        init = f(init, data[i]);
      }
    }

    return init;
  }

#ifdef VTKM_OPENMP_USE_NATIVE_REDUCTION

// Specialize for vtkm functors with OpenMP special cases:
#define VTKM_OPENMP_SPECIALIZE_REDUCE1(FunctorType, PragmaString)                                  \
  template <typename PortalT, typename ReturnType>                                                 \
  static ReturnType Execute(                                                                       \
    PortalT portal, ReturnType value, FunctorType functorIn, std::true_type)                       \
  {                                                                                                \
    const vtkm::Id numValues = portal.GetNumberOfValues();                                         \
    internal::WrappedBinaryOperator<ReturnType, FunctorType> f(functorIn);                         \
    _Pragma(#PragmaString) for (vtkm::Id i = 0; i < numValues; ++i)                                \
    {                                                                                              \
      value = f(value, portal.Get(i));                                                             \
    }                                                                                              \
    return value;                                                                                  \
  }

// Constructing the pragma string inside the _Pragma call doesn't work so
// we jump through a hoop:
#define VTKM_OPENMP_SPECIALIZE_REDUCE(FunctorType, Operator)                                       \
  VTKM_OPENMP_SPECIALIZE_REDUCE1(FunctorType, "omp parallel for reduction(" #Operator ":value)")

  // + (Add, Sum)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::Add, +)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::Sum, +)
  // * (Multiply, Product)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::Multiply, *)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::Product, *)
  // - (Subtract)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::Subtract, -)
  // & (BitwiseAnd)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::BitwiseAnd, &)
  // | (BitwiseOr)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::BitwiseOr, |)
  // ^ (BitwiseXor)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::BitwiseXor, ^)
  // && (LogicalAnd)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::LogicalAnd, &&)
  // || (LogicalOr)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::LogicalOr, ||)
  // min (Minimum)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::Minimum, min)
  // max (Maximum)
  VTKM_OPENMP_SPECIALIZE_REDUCE(vtkm::Maximum, max)

#undef VTKM_OPENMP_SPECIALIZE_REDUCE
#undef VTKM_OPENMP_SPECIALIZE_REDUCE1

#endif // VTKM_OPENMP_USE_NATIVE_REDUCTION
};

template <typename KeysInArray,
          typename ValuesInArray,
          typename KeysOutArray,
          typename ValuesOutArray,
          typename BinaryFunctor>
void ReduceByKeyHelper(KeysInArray keysInArray,
                       ValuesInArray valuesInArray,
                       KeysOutArray keysOutArray,
                       ValuesOutArray valuesOutArray,
                       BinaryFunctor functor)
{
  using KeyType = typename KeysInArray::ValueType;
  using ValueType = typename ValuesInArray::ValueType;

  const vtkm::Id numValues = keysInArray.GetNumberOfValues();
  auto keysInPortal = keysInArray.PrepareForInput(DeviceAdapterTagOpenMP());
  auto valuesInPortal = valuesInArray.PrepareForInput(DeviceAdapterTagOpenMP());
  auto keysIn = vtkm::cont::ArrayPortalToIteratorBegin(keysInPortal);
  auto valuesIn = vtkm::cont::ArrayPortalToIteratorBegin(valuesInPortal);

  auto keysOutPortal = keysOutArray.PrepareForOutput(numValues, DeviceAdapterTagOpenMP());
  auto valuesOutPortal = valuesOutArray.PrepareForOutput(numValues, DeviceAdapterTagOpenMP());
  auto keysOut = vtkm::cont::ArrayPortalToIteratorBegin(keysOutPortal);
  auto valuesOut = vtkm::cont::ArrayPortalToIteratorBegin(valuesOutPortal);

  internal::WrappedBinaryOperator<ValueType, BinaryFunctor> f(functor);
  vtkm::Id outIdx = 0;

  VTKM_OPENMP_DIRECTIVE(parallel default(none) firstprivate(keysIn, valuesIn, keysOut, valuesOut, f)
                          shared(outIdx))
  {
    int tid = omp_get_thread_num();
    int numThreads = omp_get_num_threads();

    // Determine bounds for this thread's scan operation:
    vtkm::Id chunkSize = (numValues + numThreads - 1) / numThreads;
    vtkm::Id scanIdx = std::min(tid * chunkSize, numValues);
    vtkm::Id scanEnd = std::min(scanIdx + chunkSize, numValues);

    auto threadKeysBegin = keysOut + scanIdx;
    auto threadValuesBegin = valuesOut + scanIdx;
    auto threadKey = threadKeysBegin;
    auto threadValue = threadValuesBegin;

    // Reduce each thread's partition:
    KeyType rangeKey;
    ValueType rangeValue;
    for (;;)
    {
      if (scanIdx < scanEnd)
      {
        rangeKey = keysIn[scanIdx];
        rangeValue = valuesIn[scanIdx];
        ++scanIdx;

        // Locate end of current range:
        while (scanIdx < scanEnd && static_cast<KeyType>(keysIn[scanIdx]) == rangeKey)
        {
          rangeValue = f(rangeValue, valuesIn[scanIdx]);
          ++scanIdx;
        }

        *threadKey = rangeKey;
        *threadValue = rangeValue;
        ++threadKey;
        ++threadValue;
      }
      else
      {
        break;
      }
    }

    if (tid == 0)
    {
      outIdx = static_cast<vtkm::Id>(threadKey - threadKeysBegin);
    }

    // Combine the reduction results. Skip tid == 0, since it's already in
    // the correct location:
    for (int i = 1; i < numThreads; ++i)
    {

      // This barrier ensures that:
      // 1) Threads remain synchronized through this final reduction loop.
      // 2) The outIdx variable is initialized by thread 0.
      // 3) All threads have reduced their partitions.
      VTKM_OPENMP_DIRECTIVE(barrier)

      if (tid == i)
      {
        // Check if the previous thread's last key matches our first:
        if (outIdx > 0 && threadKeysBegin < threadKey && keysOut[outIdx - 1] == *threadKeysBegin)
        {
          valuesOut[outIdx - 1] = f(valuesOut[outIdx - 1], *threadValuesBegin);
          ++threadKeysBegin;
          ++threadValuesBegin;
        }

        // Copy reduced partition to final location (if needed)
        if (threadKeysBegin < threadKey && threadKeysBegin != keysOut + outIdx)
        {
          std::copy(threadKeysBegin, threadKey, keysOut + outIdx);
          std::copy(threadValuesBegin, threadValue, valuesOut + outIdx);
        }

        outIdx += static_cast<vtkm::Id>(threadKey - threadKeysBegin);

      } // end tid == i
    }   // end combine reduction
  }     // end parallel

  keysOutArray.Shrink(outIdx);
  valuesOutArray.Shrink(outIdx);
}

template <typename IterT, typename RawPredicateT>
struct UniqueHelper
{
  using ValueType = typename std::iterator_traits<IterT>::value_type;
  using PredicateT = internal::WrappedBinaryOperator<bool, RawPredicateT>;

  struct Node
  {
    vtkm::Id2 InputRange{ -1, -1 };
    vtkm::Id2 OutputRange{ -1, -1 };

    // Pad the node out to the size of a cache line to prevent false sharing:
    static constexpr size_t DataSize = 2 * sizeof(vtkm::Id2);
    static constexpr size_t NumCacheLines = CeilDivide<size_t>(DataSize, CACHE_LINE_SIZE);
    static constexpr size_t PaddingSize = NumCacheLines * CACHE_LINE_SIZE - DataSize;
    unsigned char Padding[PaddingSize];
  };

  IterT Data;
  vtkm::Id NumValues;
  PredicateT Predicate;
  vtkm::Id LeafSize;
  std::vector<Node> Nodes;
  size_t NextNode;

  UniqueHelper(IterT iter, vtkm::Id numValues, RawPredicateT pred)
    : Data(iter)
    , NumValues(numValues)
    , Predicate(pred)
    , LeafSize(0)
    , NextNode(0)
  {
  }

  vtkm::Id Execute()
  {
    vtkm::Id outSize = 0;

    VTKM_OPENMP_DIRECTIVE(parallel default(shared))
    {
      VTKM_OPENMP_DIRECTIVE(single)
      {
        this->Prepare();

        // Kick off task-based divide-and-conquer uniquification:
        Node* rootNode = this->AllocNode();
        rootNode->InputRange = vtkm::Id2(0, this->NumValues);
        this->Uniquify(rootNode);
        outSize = rootNode->OutputRange[1] - rootNode->OutputRange[0];
      }
    }

    return outSize;
  }

private:
  void Prepare()
  {
    // Figure out how many values each thread should handle:
    int numThreads = omp_get_num_threads();
    vtkm::Id chunksPerThread = 8;
    vtkm::Id numChunks;
    ComputeChunkSize(
      this->NumValues, numThreads, chunksPerThread, sizeof(ValueType), numChunks, this->LeafSize);

    // Compute an upper-bound of the number of nodes in the tree:
568
    std::size_t numNodes = static_cast<std::size_t>(numChunks);
569
570
571
    while (numChunks > 1)
    {
      numChunks = (numChunks + 1) / 2;
572
      numNodes += static_cast<std::size_t>(numChunks);
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
    }
    this->Nodes.resize(numNodes);
    this->NextNode = 0;
  }

  Node* AllocNode()
  {
    size_t nodeIdx;

// GCC emits a false positive "value computed but not used" for this block:
#pragma GCC diagnostic push
#pragma GCC diagnostic ignored "-Wunused-value"

    VTKM_OPENMP_DIRECTIVE(atomic capture)
    {
      nodeIdx = this->NextNode;
      ++this->NextNode;
    }

#pragma GCC diagnostic pop

    VTKM_ASSERT(nodeIdx < this->Nodes.size());

    return &this->Nodes[nodeIdx];
  }

  bool IsLeaf(const vtkm::Id2& range) { return (range[1] - range[0]) <= this->LeafSize; }

  // Not an strict midpoint, but ensures that the first range will always be
  // a multiple of the leaf size.
  vtkm::Id ComputeMidpoint(const vtkm::Id2& range)
  {
    const vtkm::Id n = range[1] - range[0];
    const vtkm::Id np = this->LeafSize;

    return CeilDivide(n / 2, np) * np + range[0];
  }

  void Uniquify(Node* node)
  {
    if (!this->IsLeaf(node->InputRange))
    {
      vtkm::Id midpoint = this->ComputeMidpoint(node->InputRange);

      Node* right = this->AllocNode();
      Node* left = this->AllocNode();

      right->InputRange = vtkm::Id2(midpoint, node->InputRange[1]);

      // Intel compilers seem to have trouble following the 'this' pointer
      // when launching tasks, resulting in a corrupt task environment.
      // Explicitly copying the pointer into a local variable seems to fix this.
      auto explicitThis = this;

      VTKM_OPENMP_DIRECTIVE(taskgroup)
      {
        VTKM_OPENMP_DIRECTIVE(task) { explicitThis->Uniquify(right); }

        left->InputRange = vtkm::Id2(node->InputRange[0], midpoint);
        this->Uniquify(left);

      } // end taskgroup. Both sides of the tree will be completed here.

      // Combine the ranges in the left side:
      if (this->Predicate(this->Data[left->OutputRange[1] - 1], this->Data[right->OutputRange[0]]))
      {
        ++right->OutputRange[0];
      }

      vtkm::Id numVals = right->OutputRange[1] - right->OutputRange[0];
      DoCopy(this->Data + right->OutputRange[0], this->Data + left->OutputRange[1], numVals);

      node->OutputRange[0] = left->OutputRange[0];
      node->OutputRange[1] = left->OutputRange[1] + numVals;
    }
    else
    {
      auto start = this->Data + node->InputRange[0];
      auto end = this->Data + node->InputRange[1];
      end = std::unique(start, end, this->Predicate);
      node->OutputRange[0] = node->InputRange[0];
      node->OutputRange[1] = node->InputRange[0] + static_cast<vtkm::Id>(end - start);
    }
  }
};
}
}
} // end namespace vtkm::cont::openmp

#endif // vtk_m_cont_openmp_internal_FunctorsOpenMP_h