ITK  5.4.0
Insight Toolkit
itkMultiThreaderBase.h
Go to the documentation of this file.
1/*=========================================================================
2 *
3 * Copyright NumFOCUS
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License");
6 * you may not use this file except in compliance with the License.
7 * You may obtain a copy of the License at
8 *
9 * https://www.apache.org/licenses/LICENSE-2.0.txt
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 *
17 *=========================================================================*/
18/*=========================================================================
19 *
20 * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21 *
22 * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23 *
24 * For complete copyright, license and disclaimer of warranty information
25 * please refer to the NOTICE file at the top of the ITK source tree.
26 *
27 *=========================================================================*/
28#ifndef itkMultiThreaderBase_h
29#define itkMultiThreaderBase_h
30
31#include "itkObject.h"
32#include "itkThreadSupport.h"
33#include "itkObjectFactory.h"
34#include "itkIntTypes.h"
35#include "itkImageRegion.h"
36#include "itkImageIORegion.h"
37#include "itkSingletonMacro.h"
38#include <atomic>
39#include <functional>
40#include <thread>
41#include "itkProgressReporter.h"
42
43
44namespace itk
45{
46
54{
55public:
61 enum class Threader : int8_t
62 {
63 Platform = 0,
64 First = Platform,
65 Pool,
66 TBB,
67 Last = TBB,
68 Unknown = -1
69 };
70
74 enum class ThreadExitCode : uint8_t
75 {
76 SUCCESS,
77 ITK_EXCEPTION,
78 ITK_PROCESS_ABORTED_EXCEPTION,
79 STD_EXCEPTION,
80 UNKNOWN
81 };
82};
83// Define how to print enumeration
84extern ITKCommon_EXPORT std::ostream &
85 operator<<(std::ostream & out, const MultiThreaderBaseEnums::Threader value);
86extern ITKCommon_EXPORT std::ostream &
87 operator<<(std::ostream & out, const MultiThreaderBaseEnums::ThreadExitCode value);
88
103struct MultiThreaderBaseGlobals;
104class ProcessObject;
105
106class ITKCommon_EXPORT MultiThreaderBase : public Object
107{
108public:
109 ITK_DISALLOW_COPY_AND_MOVE(MultiThreaderBase);
110
116
118 static Pointer
120
122 itkOverrideGetNameOfClassMacro(MultiThreaderBase);
123
127 virtual void
129 itkGetConstMacro(MaximumNumberOfThreads, ThreadIdType);
135 virtual void
137 itkGetConstMacro(NumberOfWorkUnits, ThreadIdType);
140 virtual void
141 SetUpdateProgress(bool updates);
142 itkGetConstMacro(UpdateProgress, bool);
143
149 static void
151 static ThreadIdType
160 itkLegacyMacro(static void SetGlobalDefaultUseThreadPool(const bool GlobalDefaultUseThreadPool));
161 itkLegacyMacro(static bool GetGlobalDefaultUseThreadPool());
165#if !defined(ITK_LEGACY_REMOVE)
166 using ThreaderType = ThreaderEnum;
167 static constexpr ThreaderEnum Platform = ThreaderEnum::Platform;
168 static constexpr ThreaderEnum First = ThreaderEnum::First;
169 static constexpr ThreaderEnum Pool = ThreaderEnum::Pool;
170 static constexpr ThreaderEnum TBB = ThreaderEnum::TBB;
171 static constexpr ThreaderEnum Last = ThreaderEnum::Last;
172 static constexpr ThreaderEnum Unknown = ThreaderEnum::Unknown;
173#endif
174
176 static ThreaderEnum
177 ThreaderTypeFromString(std::string threaderString);
178
180 static std::string
182 {
183 switch (threader)
184 {
185 case ThreaderEnum::Platform:
186 return "Platform";
187 case ThreaderEnum::Pool:
188 return "Pool";
189 case ThreaderEnum::TBB:
190 return "TBB";
191 case ThreaderEnum::Unknown:
192 default:
193 return "Unknown";
194 }
195 }
209 static void
211 static ThreaderEnum
219 static void
221 static ThreadIdType
225#if !defined(ITK_LEGACY_REMOVE)
228 itkLegacyMacro(virtual void SetNumberOfThreads(ThreadIdType numberOfThreads))
229 {
230 this->SetMaximumNumberOfThreads(numberOfThreads);
231 this->SetNumberOfWorkUnits(this->GetMaximumNumberOfThreads()); // Might be clamped
232 }
233 itkLegacyMacro(virtual ThreadIdType GetNumberOfThreads()) { return this->GetNumberOfWorkUnits(); }
244 // clang-format off
245ITK_GCC_PRAGMA_DIAG_PUSH()
246ITK_GCC_PRAGMA_DIAG(ignored "-Wattributes")
247INTEL_PRAGMA_WARN_PUSH
248INTEL_SUPPRESS_warning_1292
249 // clang-format on
250# ifdef ITK_LEGACY_SILENT
251 struct ThreadInfoStruct
252# else
253 struct [[deprecated("Use WorkUnitInfo, ThreadInfoStruct is deprecated since ITK 5.0")]] ThreadInfoStruct
254# endif
255 // clang-format off
256INTEL_PRAGMA_WARN_POP
257 // clang-format on
258 {
259 ThreadIdType ThreadID;
260 ThreadIdType NumberOfThreads;
261 void * UserData;
262 ThreadFunctionType ThreadFunction;
263 enum
264 {
265 SUCCESS,
266 ITK_EXCEPTION,
267 ITK_PROCESS_ABORTED_EXCEPTION,
268 STD_EXCEPTION,
269 UNKNOWN
271 };
272 // clang-format off
273ITK_GCC_PRAGMA_DIAG_POP()
274 // clang-format on
275#endif // ITK_LEGACY_REMOVE
284 struct WorkUnitInfo
285 {
288 void * UserData;
292#if !defined(ITK_LEGACY_REMOVE)
293 using ThreadExitCodeType = ThreadExitCodeEnum;
294 static constexpr ThreadExitCodeEnum SUCCESS = ThreadExitCodeEnum::SUCCESS;
295 static constexpr ThreadExitCodeEnum ITK_EXCEPTION = ThreadExitCodeEnum::ITK_EXCEPTION;
296 static constexpr ThreadExitCodeEnum ITK_PROCESS_ABORTED_EXCEPTION =
297 ThreadExitCodeEnum::ITK_PROCESS_ABORTED_EXCEPTION;
298 static constexpr ThreadExitCodeEnum STD_EXCEPTION = ThreadExitCodeEnum::STD_EXCEPTION;
299 static constexpr ThreadExitCodeEnum UNKNOWN = ThreadExitCodeEnum::UNKNOWN;
300#endif
301 };
308 virtual void
310
315 virtual void
317
319 void
321
322#ifndef ITK_FUTURE_LEGACY_REMOVE
323 // `TemplatedThreadingFunctorType` was previously used to declare the `funcP` parameter of `ParallelizeImageRegion`
324 // and `ParallelizeImageRegionRestrictDirection` template member functions.
325 template <unsigned int VDimension>
326 using TemplatedThreadingFunctorType [[deprecated]] = std::function<void(const ImageRegion<VDimension> &)>;
327#endif
328
329 using ThreadingFunctorType = std::function<void(const IndexValueType index[], const SizeValueType size[])>;
330 using ArrayThreadingFunctorType = std::function<void(SizeValueType)>;
331
337 virtual void
339 SizeValueType lastIndexPlus1,
341 ProcessObject * filter);
342
348 template <unsigned int VDimension, typename TFunction>
349 ITK_TEMPLATE_EXPORT void
350 ParallelizeImageRegion(const ImageRegion<VDimension> & requestedRegion, TFunction funcP, ProcessObject * filter)
351 {
352 this->ParallelizeImageRegion(
353 VDimension,
354 requestedRegion.GetIndex().m_InternalArray,
355 requestedRegion.GetSize().m_InternalArray,
356 [&funcP](const IndexValueType index[], const SizeValueType size[]) {
357 ImageRegion<VDimension> region;
358 for (unsigned int d = 0; d < VDimension; ++d)
359 {
360 region.SetIndex(d, index[d]);
361 region.SetSize(d, size[d]);
362 }
363 funcP(region);
364 },
365 filter);
366 }
372 template <unsigned int VDimension, typename TFunction>
373 ITK_TEMPLATE_EXPORT void
374 ParallelizeImageRegionRestrictDirection(unsigned int restrictedDirection,
375 const ImageRegion<VDimension> & requestedRegion,
376 TFunction funcP,
377 ProcessObject * filter)
378 {
379 if constexpr (VDimension <= 1) // Cannot split, no parallelization
380 {
381
382 ProgressReporter progress(filter, 0, requestedRegion.GetNumberOfPixels());
383 funcP(requestedRegion);
384 }
385 else // Can split, parallelize!
386 {
387 constexpr unsigned int SplitDimension = VDimension - 1;
388
389 Index<SplitDimension> splitIndex{};
390 Size<SplitDimension> splitSize{};
391 for (unsigned int splitDimension = 0, dimension = 0; dimension < VDimension; ++dimension)
392 {
393 if (dimension != restrictedDirection)
394 {
395 splitIndex[splitDimension] = requestedRegion.GetIndex(dimension);
396 splitSize[splitDimension] = requestedRegion.GetSize(dimension);
397 ++splitDimension;
398 }
399 }
400
401 this->ParallelizeImageRegion(
402 SplitDimension,
403 splitIndex.m_InternalArray,
404 splitSize.m_InternalArray,
405 [restrictedDirection, &requestedRegion, &funcP](const IndexValueType index[], const SizeValueType size[]) {
406 ImageRegion<VDimension> restrictedRequestedRegion;
407 restrictedRequestedRegion.SetIndex(restrictedDirection, requestedRegion.GetIndex(restrictedDirection));
408 restrictedRequestedRegion.SetSize(restrictedDirection, requestedRegion.GetSize(restrictedDirection));
409 for (unsigned int splitDimension = 0, dimension = 0; dimension < VDimension; ++dimension)
410 {
411 if (dimension != restrictedDirection)
412 {
413 restrictedRequestedRegion.SetIndex(dimension, index[splitDimension]);
414 restrictedRequestedRegion.SetSize(dimension, size[splitDimension]);
415 ++splitDimension;
416 }
417 }
418 funcP(restrictedRequestedRegion);
419 },
420 filter);
421 }
422 }
423
426 virtual void
427 ParallelizeImageRegion(unsigned int dimension,
428 const IndexValueType index[],
429 const SizeValueType size[],
431 ProcessObject * filter);
432
433protected:
436 void
437 PrintSelf(std::ostream & os, Indent indent) const override;
438
440 {
445 };
446
449
451 {
453 unsigned int dimension;
457 };
458
461
463 ThreadIdType m_NumberOfWorkUnits{};
464
473 ThreadIdType m_MaximumNumberOfThreads{};
474
482 SingleMethodProxy(void * arg);
483
485 ThreadFunctionType m_SingleMethod{ nullptr };
486
488 void * m_SingleData{ nullptr };
489
490private:
491 static void
493 static ThreaderEnum
495
497 itkGetGlobalDeclarationMacro(MultiThreaderBaseGlobals, PimplGlobals);
498
499 std::atomic<bool> m_UpdateProgress{ true };
500
501 static MultiThreaderBaseGlobals * m_PimplGlobals;
505 friend class ProcessObject;
506
508 static ThreadIdType
510};
511
512
513} // namespace itk
514
515#endif
SizeValueType GetNumberOfPixels() const
const IndexType & GetIndex() const
const SizeType & GetSize() const
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
enums for MultiThreaderBase
A class for performing multithreaded execution.
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ParallelizeArrayHelper(void *arg)
ITK_TEMPLATE_EXPORT void ParallelizeImageRegion(const ImageRegion< VDimension > &requestedRegion, TFunction funcP, ProcessObject *filter)
ITK_TEMPLATE_EXPORT void ParallelizeImageRegionRestrictDirection(unsigned int restrictedDirection, const ImageRegion< VDimension > &requestedRegion, TFunction funcP, ProcessObject *filter)
virtual void SetMaximumNumberOfThreads(ThreadIdType numberOfThreads)
void SetSingleMethodAndExecute(ThreadFunctionType func, void *data)
static std::string ThreaderTypeToString(ThreaderEnum threader)
static ThreaderEnum GetGlobalDefaultThreaderPrivate()
static Pointer New()
static void SetGlobalDefaultThreaderPrivate(ThreaderEnum threaderType)
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION SingleMethodProxy(void *arg)
static ThreaderEnum GetGlobalDefaultThreader()
static void SetGlobalDefaultThreader(ThreaderEnum threaderType)
static MultiThreaderBaseGlobals * m_PimplGlobals
static ThreaderEnum ThreaderTypeFromString(std::string threaderString)
static ThreadIdType GetGlobalMaximumNumberOfThreads()
virtual void SetNumberOfWorkUnits(ThreadIdType numberOfWorkUnits)
virtual void SetUpdateProgress(bool updates)
static ThreadIdType GetGlobalDefaultNumberOfThreadsByPlatform()
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void SetSingleMethod(ThreadFunctionType, void *data)=0
itkGetGlobalDeclarationMacro(MultiThreaderBaseGlobals, PimplGlobals)
virtual void ParallelizeImageRegion(unsigned int dimension, const IndexValueType index[], const SizeValueType size[], ThreadingFunctorType funcP, ProcessObject *filter)
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ParallelizeImageRegionHelper(void *arg)
static ThreadIdType GetGlobalDefaultNumberOfThreads()
itkLegacyMacro(static void SetGlobalDefaultUseThreadPool(const bool GlobalDefaultUseThreadPool))
itkLegacyMacro(static bool GetGlobalDefaultUseThreadPool())
virtual void SingleMethodExecute()=0
std::function< void(SizeValueType)> ArrayThreadingFunctorType
~MultiThreaderBase() override
std::function< void(const IndexValueType index[], const SizeValueType size[])> ThreadingFunctorType
static void SetGlobalMaximumNumberOfThreads(ThreadIdType val)
static void SetGlobalDefaultNumberOfThreads(ThreadIdType val)
virtual void ParallelizeArray(SizeValueType firstIndex, SizeValueType lastIndexPlus1, ArrayThreadingFunctorType aFunc, ProcessObject *filter)
Base class for most ITK classes.
Definition: itkObject.h:62
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Implements progress tracking for a filter.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
long IndexValueType
Definition: itkIntTypes.h:90
void(*)(void *) ThreadFunctionType
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:216
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::ITK_THREAD_RETURN_TYPE ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:71
IndexValueType m_InternalArray[VDimension]
Definition: itkIndex.h:293
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:72
SizeValueType m_InternalArray[VDimension]
Definition: itkSize.h:245