forked from Kitware/VTK
-
Notifications
You must be signed in to change notification settings - Fork 0
/
vtkSMPTools.h
602 lines (551 loc) · 18.5 KB
/
vtkSMPTools.h
1
2
3
4
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
33
34
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
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
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
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
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
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
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
568
569
570
571
572
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
/*=========================================================================
Program: Visualization Toolkit
Module: vtkSMPTools.h
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm 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.
=========================================================================*/
/**
* @class vtkSMPTools
* @brief A set of parallel (multi-threaded) utility functions.
*
* vtkSMPTools provides a set of utility functions that can
* be used to parallelize parts of VTK code using multiple threads.
* There are several back-end implementations of parallel functionality
* (currently Sequential, TBB, OpenMP and STDThread) that actual execution is
* delegated to.
*
* @sa
* vtkSMPThreadLocal
* vtkSMPThreadLocalObject
*/
#ifndef vtkSMPTools_h
#define vtkSMPTools_h
#include "vtkCommonCoreModule.h" // For export macro
#include "vtkObject.h"
#include "SMP/Common/vtkSMPToolsAPI.h"
#include "vtkSMPThreadLocal.h" // For Initialized
#include <functional> // For std::function
#include <iterator> // For std::iterator
#include <type_traits> // For std:::enable_if
#ifndef DOXYGEN_SHOULD_SKIP_THIS
#ifndef __VTK_WRAP__
namespace vtk
{
namespace detail
{
namespace smp
{
template <typename T>
class vtkSMPTools_Has_Initialize
{
typedef char (&no_type)[1];
typedef char (&yes_type)[2];
template <typename U, void (U::*)()>
struct V
{
};
template <typename U>
static yes_type check(V<U, &U::Initialize>*);
template <typename U>
static no_type check(...);
public:
static bool const value = sizeof(check<T>(nullptr)) == sizeof(yes_type);
};
template <typename T>
class vtkSMPTools_Has_Initialize_const
{
typedef char (&no_type)[1];
typedef char (&yes_type)[2];
template <typename U, void (U::*)() const>
struct V
{
};
template <typename U>
static yes_type check(V<U, &U::Initialize>*);
template <typename U>
static no_type check(...);
public:
static bool const value = sizeof(check<T>(0)) == sizeof(yes_type);
};
template <typename Functor, bool Init>
struct vtkSMPTools_FunctorInternal;
template <typename Functor>
struct vtkSMPTools_FunctorInternal<Functor, false>
{
Functor& F;
vtkSMPTools_FunctorInternal(Functor& f)
: F(f)
{
}
void Execute(vtkIdType first, vtkIdType last) { this->F(first, last); }
void For(vtkIdType first, vtkIdType last, vtkIdType grain)
{
auto& SMPToolsAPI = vtkSMPToolsAPI::GetInstance();
SMPToolsAPI.For(first, last, grain, *this);
}
vtkSMPTools_FunctorInternal<Functor, false>& operator=(
const vtkSMPTools_FunctorInternal<Functor, false>&);
vtkSMPTools_FunctorInternal<Functor, false>(const vtkSMPTools_FunctorInternal<Functor, false>&);
};
template <typename Functor>
struct vtkSMPTools_FunctorInternal<Functor, true>
{
Functor& F;
vtkSMPThreadLocal<unsigned char> Initialized;
vtkSMPTools_FunctorInternal(Functor& f)
: F(f)
, Initialized(0)
{
}
void Execute(vtkIdType first, vtkIdType last)
{
unsigned char& inited = this->Initialized.Local();
if (!inited)
{
this->F.Initialize();
inited = 1;
}
this->F(first, last);
}
void For(vtkIdType first, vtkIdType last, vtkIdType grain)
{
auto& SMPToolsAPI = vtkSMPToolsAPI::GetInstance();
SMPToolsAPI.For(first, last, grain, *this);
this->F.Reduce();
}
vtkSMPTools_FunctorInternal<Functor, true>& operator=(
const vtkSMPTools_FunctorInternal<Functor, true>&);
vtkSMPTools_FunctorInternal<Functor, true>(const vtkSMPTools_FunctorInternal<Functor, true>&);
};
template <typename Functor>
class vtkSMPTools_Lookup_For
{
static bool const init = vtkSMPTools_Has_Initialize<Functor>::value;
public:
typedef vtkSMPTools_FunctorInternal<Functor, init> type;
};
template <typename Functor>
class vtkSMPTools_Lookup_For<Functor const>
{
static bool const init = vtkSMPTools_Has_Initialize_const<Functor>::value;
public:
typedef vtkSMPTools_FunctorInternal<Functor const, init> type;
};
template <typename Iterator, typename Functor, bool Init>
struct vtkSMPTools_RangeFunctor;
template <typename Iterator, typename Functor>
struct vtkSMPTools_RangeFunctor<Iterator, Functor, false>
{
Functor& F;
Iterator& Begin;
vtkSMPTools_RangeFunctor(Iterator& begin, Functor& f)
: F(f)
, Begin(begin)
{
}
void operator()(vtkIdType first, vtkIdType last)
{
Iterator itFirst(Begin);
std::advance(itFirst, first);
Iterator itLast(itFirst);
std::advance(itLast, last - first);
this->F(itFirst, itLast);
}
};
template <typename Iterator, typename Functor>
struct vtkSMPTools_RangeFunctor<Iterator, Functor, true>
{
Functor& F;
Iterator& Begin;
vtkSMPTools_RangeFunctor(Iterator& begin, Functor& f)
: F(f)
, Begin(begin)
{
}
void Initialize() { this->F.Initialize(); }
void operator()(vtkIdType first, vtkIdType last)
{
Iterator itFirst(Begin);
std::advance(itFirst, first);
Iterator itLast(itFirst);
std::advance(itLast, last - first);
this->F(itFirst, itLast);
}
void Reduce() { this->F.Reduce(); }
};
template <typename Iterator, typename Functor>
class vtkSMPTools_Lookup_RangeFor
{
static bool const init = vtkSMPTools_Has_Initialize<Functor>::value;
public:
typedef vtkSMPTools_RangeFunctor<Iterator, Functor, init> type;
};
template <typename Iterator, typename Functor>
class vtkSMPTools_Lookup_RangeFor<Iterator, Functor const>
{
static bool const init = vtkSMPTools_Has_Initialize_const<Functor>::value;
public:
typedef vtkSMPTools_RangeFunctor<Iterator, Functor const, init> type;
};
template <typename T>
using resolvedNotInt = typename std::enable_if<!std::is_integral<T>::value, void>::type;
} // namespace smp
} // namespace detail
} // namespace vtk
#endif // __VTK_WRAP__
#endif // DOXYGEN_SHOULD_SKIP_THIS
class VTKCOMMONCORE_EXPORT vtkSMPTools
{
public:
///@{
/**
* Execute a for operation in parallel. First and last
* define the range over which to operate (which is defined
* by the operator). The operation executed is defined by
* operator() of the functor object. The grain gives the parallel
* engine a hint about the coarseness over which to parallelize
* the function (as defined by last-first of each execution of
* operator() ).
*/
template <typename Functor>
static void For(vtkIdType first, vtkIdType last, vtkIdType grain, Functor& f)
{
typename vtk::detail::smp::vtkSMPTools_Lookup_For<Functor>::type fi(f);
fi.For(first, last, grain);
}
template <typename Functor>
static void For(vtkIdType first, vtkIdType last, vtkIdType grain, Functor const& f)
{
typename vtk::detail::smp::vtkSMPTools_Lookup_For<Functor const>::type fi(f);
fi.For(first, last, grain);
}
///@}
///@{
/**
* Execute a for operation in parallel. First and last
* define the range over which to operate (which is defined
* by the operator). The operation executed is defined by
* operator() of the functor object. The grain gives the parallel
* engine a hint about the coarseness over which to parallelize
* the function (as defined by last-first of each execution of
* operator() ). Uses a default value for the grain.
*/
template <typename Functor>
static void For(vtkIdType first, vtkIdType last, Functor& f)
{
vtkSMPTools::For(first, last, 0, f);
}
template <typename Functor>
static void For(vtkIdType first, vtkIdType last, Functor const& f)
{
vtkSMPTools::For(first, last, 0, f);
}
///@}
///@{
/**
* Execute a for operation in parallel. Begin and end iterators
* define the range over which to operate (which is defined
* by the operator). The operation executed is defined by
* operator() of the functor object. The grain gives the parallel
* engine a hint about the coarseness over which to parallelize
* the function (as defined by last-first of each execution of
* operator() ).
*
* Usage example:
* \code
* template<class IteratorT>
* class ExampleFunctor
* {
* void operator()(IteratorT begin, IteratorT end)
* {
* for (IteratorT it = begin; it != end; ++it)
* {
* // Do stuff
* }
* }
* };
* ExampleFunctor<std::set<int>::iterator> worker;
* vtkSMPTools::For(container.begin(), container.end(), 5, worker);
* \endcode
*
* Lambda are also supported through Functor const&
* function overload:
* \code
* vtkSMPTools::For(container.begin(), container.end(), 5,
* [](std::set<int>::iterator begin, std::set<int>::iterator end) {
* // Do stuff
* });
* \endcode
*/
template <typename Iter, typename Functor>
static vtk::detail::smp::resolvedNotInt<Iter> For(
Iter begin, Iter end, vtkIdType grain, Functor& f)
{
vtkIdType size = std::distance(begin, end);
typename vtk::detail::smp::vtkSMPTools_Lookup_RangeFor<Iter, Functor>::type fi(begin, f);
vtkSMPTools::For(0, size, grain, fi);
}
template <typename Iter, typename Functor>
static vtk::detail::smp::resolvedNotInt<Iter> For(
Iter begin, Iter end, vtkIdType grain, Functor const& f)
{
vtkIdType size = std::distance(begin, end);
typename vtk::detail::smp::vtkSMPTools_Lookup_RangeFor<Iter, Functor const>::type fi(begin, f);
vtkSMPTools::For(0, size, grain, fi);
}
///@}
///@{
/**
* Execute a for operation in parallel. Begin and end iterators
* define the range over which to operate (which is defined
* by the operator). The operation executed is defined by
* operator() of the functor object. Uses a default value
* for the grain.
*
* Usage example:
* \code
* template<class IteratorT>
* class ExampleFunctor
* {
* void operator()(IteratorT begin, IteratorT end)
* {
* for (IteratorT it = begin; it != end; ++it)
* {
* // Do stuff
* }
* }
* };
* ExampleFunctor<std::set<int>::iterator> worker;
* vtkSMPTools::For(container.begin(), container.end(), worker);
* \endcode
*
* Lambda are also supported through Functor const&
* function overload:
* \code
* vtkSMPTools::For(container.begin(), container.end(),
* [](std::set<int>::iterator begin, std::set<int>::iterator end) {
* // Do stuff
* });
* \endcode
*/
template <typename Iter, typename Functor>
static vtk::detail::smp::resolvedNotInt<Iter> For(Iter begin, Iter end, Functor& f)
{
vtkSMPTools::For(begin, end, 0, f);
}
template <typename Iter, typename Functor>
static vtk::detail::smp::resolvedNotInt<Iter> For(Iter begin, Iter end, Functor const& f)
{
vtkSMPTools::For(begin, end, 0, f);
}
///@}
/**
* Get the backend in use.
*/
static const char* GetBackend();
/**
* /!\ This method is not thread safe.
* Change the backend in use.
* The options can be: "Sequential", "STDThread", "TBB" or "OpenMP"
*
* VTK_SMP_BACKEND_IN_USE env variable can also be used to set the default SMPTools
* backend, in that case SetBackend() doesn't need to be called.
* The backend selected with SetBackend() have the priority over VTK_SMP_BACKEND_IN_USE.
*
* SetBackend() will return true if the backend was found and available.
*/
static bool SetBackend(const char* backend);
/**
* /!\ This method is not thread safe.
* Initialize the underlying libraries for execution. This is
* not required as it is automatically defined by the libaries.
* However, it can be used to control the maximum number of thread used.
* Make sure to call it before the parallel operation.
*
* If Initialize is called without argument it will reset
* to the maximum number of threads or use the VTK_SMP_MAX_THREADS
* env variable if it is defined.
*
* Note: If VTK_SMP_MAX_THREADS env variable is defined the SMPTools will try
* to use it to set the maximum number of threads. Initialize() doesn't
* need to be called.
*/
static void Initialize(int numThreads = 0);
/**
* Get the estimated number of threads being used by the backend.
* This should be used as just an estimate since the number of threads may
* vary dynamically and a particular task may not be executed on all the
* available threads.
*/
static int GetEstimatedNumberOfThreads();
/**
* /!\ This method is not thread safe.
* If true enable nested parallelism for underlying backends.
* When enabled the comportement is different for each backend:
* - TBB support nested parallelism by default.
* - For OpenMP, we set `omp_set_nested` to true so that it is supported.
* - STDThread support nested parallelism by creating new threads pools.
* - For Sequential nothing changes.
*
* Default to true
*/
static void SetNestedParallelism(bool isNested);
/**
* Get true if the nested parallelism is enabled.
*/
static bool GetNestedParallelism();
/**
* Return true if it is called from a parallel scope.
*/
static bool IsParallelScope();
/**
* Structure used to specify configuration for LocalScope() method.
* Several parameters can be configured:
* - MaxNumberOfThreads set the maximum number of threads.
* - Backend set a specific SMPTools backend.
* - NestedParallelism, if true enable nested parallelism.
*/
struct Config
{
int MaxNumberOfThreads = 0;
std::string Backend = vtk::detail::smp::vtkSMPToolsAPI::GetInstance().GetBackend();
bool NestedParallelism = true;
Config() {}
Config(int maxNumberOfThreads)
: MaxNumberOfThreads(maxNumberOfThreads)
{
}
Config(std::string backend)
: Backend(backend)
{
}
Config(bool nestedParallelism)
: NestedParallelism(nestedParallelism)
{
}
Config(int maxNumberOfThreads, std::string backend, bool nestedParallelism)
: MaxNumberOfThreads(maxNumberOfThreads)
, Backend(backend)
, NestedParallelism(nestedParallelism)
{
}
#ifndef DOXYGEN_SHOULD_SKIP_THIS
Config(vtk::detail::smp::vtkSMPToolsAPI& API)
: MaxNumberOfThreads(API.GetInternalDesiredNumberOfThread())
, Backend(API.GetBackend())
, NestedParallelism(API.GetNestedParallelism())
{
}
#endif // DOXYGEN_SHOULD_SKIP_THIS
};
/**
* /!\ This method is not thread safe.
* Change the number of threads locally within this scope and call a functor which
* should contains a vtkSMPTools method.
*
* Usage example:
* \code
* vtkSMPTools::LocalScope(
* vtkSMPTools::Config{ 4, "OpenMP", false }, [&]() { vtkSMPTools::For(0, size, worker); });
* \endcode
*/
template <typename T>
static void LocalScope(Config const& config, T&& lambda)
{
auto& SMPToolsAPI = vtk::detail::smp::vtkSMPToolsAPI::GetInstance();
SMPToolsAPI.LocalScope<vtkSMPTools::Config>(config, lambda);
}
/**
* A convenience method for transforming data. It is a drop in replacement for
* std::transform(), it does a unary operation on the input ranges. The data array must have the
* same length. The performed transformation is defined by operator() of the functor object.
*
* Usage example with vtkDataArray:
* \code
* const auto range0 = vtk::DataArrayValueRange<1>(array0);
* auto range1 = vtk::DataArrayValueRange<1>(array1);
* vtkSMPTools::Transform(
* range0.cbegin(), range0.cend(), range1.begin(), [](double x) { return x - 1; });
* \endcode
*
* Please visit vtkDataArrayRange.h documentation for more information and optimisation.
*/
template <typename InputIt, typename OutputIt, typename Functor>
static void Transform(InputIt inBegin, InputIt inEnd, OutputIt outBegin, Functor transform)
{
auto& SMPToolsAPI = vtk::detail::smp::vtkSMPToolsAPI::GetInstance();
SMPToolsAPI.Transform(inBegin, inEnd, outBegin, transform);
}
/**
* A convenience method for transforming data. It is a drop in replacement for
* std::transform(), it does a binary operation on the input ranges. The data array must have the
* same length. The performed transformation is defined by operator() of the functor object.
*
* Usage example with vtkDataArray:
* \code
* const auto range0 = vtk::DataArrayValueRange<1>(array0);
* auto range1 = vtk::DataArrayValueRange<1>(array1);
* vtkSMPTools::Transform(
* range0.cbegin(), range0.cend(), range1.cbegin(), range1.begin(),
* [](double x, double y) { return x * y; });
* \endcode
*
* Please visit vtkDataArrayRange.h documentation for more information and optimisation.
*/
template <typename InputIt1, typename InputIt2, typename OutputIt, typename Functor>
static void Transform(
InputIt1 inBegin1, InputIt1 inEnd, InputIt2 inBegin2, OutputIt outBegin, Functor transform)
{
auto& SMPToolsAPI = vtk::detail::smp::vtkSMPToolsAPI::GetInstance();
SMPToolsAPI.Transform(inBegin1, inEnd, inBegin2, outBegin, transform);
}
/**
* A convenience method for filling data. It is a drop in replacement for std::fill(),
* it assign the given value to the element in ranges.
*
* Usage example with vtkDataArray:
* \code
* // Fill range with its first tuple value
* auto range = vtk::DataArrayTupleRange<1>(array);
* const auto value = *range.begin();
* vtkSMPTools::Fill(range.begin(), range.end(), value);
* \endcode
*
* Please visit vtkDataArrayRange.h documentation for more information and optimisation.
*/
template <typename Iterator, typename T>
static void Fill(Iterator begin, Iterator end, const T& value)
{
auto& SMPToolsAPI = vtk::detail::smp::vtkSMPToolsAPI::GetInstance();
SMPToolsAPI.Fill(begin, end, value);
}
/**
* A convenience method for sorting data. It is a drop in replacement for
* std::sort(). Under the hood different methods are used. For example,
* tbb::parallel_sort is used in TBB.
*/
template <typename RandomAccessIterator>
static void Sort(RandomAccessIterator begin, RandomAccessIterator end)
{
auto& SMPToolsAPI = vtk::detail::smp::vtkSMPToolsAPI::GetInstance();
SMPToolsAPI.Sort(begin, end);
}
/**
* A convenience method for sorting data. It is a drop in replacement for
* std::sort(). Under the hood different methods are used. For example,
* tbb::parallel_sort is used in TBB. This version of Sort() takes a
* comparison class.
*/
template <typename RandomAccessIterator, typename Compare>
static void Sort(RandomAccessIterator begin, RandomAccessIterator end, Compare comp)
{
auto& SMPToolsAPI = vtk::detail::smp::vtkSMPToolsAPI::GetInstance();
SMPToolsAPI.Sort(begin, end, comp);
}
};
#endif
// VTK-HeaderTest-Exclude: vtkSMPTools.h