ITK  6.0.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);
104struct MultiThreaderBaseGlobals;
105class ProcessObject;
106
107class ITKCommon_EXPORT MultiThreaderBase : public Object
108{
109public:
110 ITK_DISALLOW_COPY_AND_MOVE(MultiThreaderBase);
111
117
119 static Pointer
121
123 itkOverrideGetNameOfClassMacro(MultiThreaderBase);
124
128 virtual void
130 itkGetConstMacro(MaximumNumberOfThreads, ThreadIdType);
136 virtual void
138 itkGetConstMacro(NumberOfWorkUnits, ThreadIdType);
141 virtual void
142 SetUpdateProgress(bool updates);
143 itkGetConstMacro(UpdateProgress, bool);
144
150 static void
152 static ThreadIdType
161 itkLegacyMacro(static void SetGlobalDefaultUseThreadPool(const bool GlobalDefaultUseThreadPool);)
162 itkLegacyMacro(static bool GetGlobalDefaultUseThreadPool();)
166#if !defined(ITK_LEGACY_REMOVE)
167 using ThreaderType = ThreaderEnum;
168 static constexpr ThreaderEnum Platform = ThreaderEnum::Platform;
169 static constexpr ThreaderEnum First = ThreaderEnum::First;
170 static constexpr ThreaderEnum Pool = ThreaderEnum::Pool;
171 static constexpr ThreaderEnum TBB = ThreaderEnum::TBB;
172 static constexpr ThreaderEnum Last = ThreaderEnum::Last;
173 static constexpr ThreaderEnum Unknown = ThreaderEnum::Unknown;
174#endif
175
177 static ThreaderEnum
178 ThreaderTypeFromString(std::string threaderString);
179
181 static std::string
183 {
184 switch (threader)
185 {
186 case ThreaderEnum::Platform:
187 return "Platform";
188 case ThreaderEnum::Pool:
189 return "Pool";
190 case ThreaderEnum::TBB:
191 return "TBB";
192 case ThreaderEnum::Unknown:
193 default:
194 return "Unknown";
195 }
196 }
210 static void
212 static ThreaderEnum
220 static void
222 static ThreadIdType
226#if !defined(ITK_LEGACY_REMOVE)
229 itkLegacyMacro(virtual void SetNumberOfThreads(ThreadIdType numberOfThreads))
230 {
231 this->SetMaximumNumberOfThreads(numberOfThreads);
232 this->SetNumberOfWorkUnits(this->GetMaximumNumberOfThreads()); // Might be clamped
233 }
234 itkLegacyMacro(virtual ThreadIdType GetNumberOfThreads())
235 {
236 return this->GetNumberOfWorkUnits();
237 }
248 // clang-format off
249ITK_GCC_PRAGMA_DIAG_PUSH()
250ITK_GCC_PRAGMA_DIAG(ignored "-Wattributes")
251INTEL_PRAGMA_WARN_PUSH
252INTEL_SUPPRESS_warning_1292
253 // clang-format on
254# ifdef ITK_LEGACY_SILENT
255 struct ThreadInfoStruct
256# else
257 struct [[deprecated("Use WorkUnitInfo, ThreadInfoStruct is deprecated since ITK 5.0")]] ThreadInfoStruct
258# endif
259 // clang-format off
260INTEL_PRAGMA_WARN_POP
261 // clang-format on
262 {
263 ThreadIdType ThreadID;
264 ThreadIdType NumberOfThreads;
265 void * UserData;
266 ThreadFunctionType ThreadFunction;
267 enum
268 {
269 SUCCESS,
270 ITK_EXCEPTION,
271 ITK_PROCESS_ABORTED_EXCEPTION,
272 STD_EXCEPTION,
273 UNKNOWN
274 } ThreadExitCode;
275 };
276 // clang-format off
277ITK_GCC_PRAGMA_DIAG_POP()
278 // clang-format on
279#endif // ITK_LEGACY_REMOVE
288 struct WorkUnitInfo
289 {
292 void * UserData;
296#if !defined(ITK_LEGACY_REMOVE)
297 using ThreadExitCodeType = ThreadExitCodeEnum;
298 static constexpr ThreadExitCodeEnum SUCCESS = ThreadExitCodeEnum::SUCCESS;
299 static constexpr ThreadExitCodeEnum ITK_EXCEPTION = ThreadExitCodeEnum::ITK_EXCEPTION;
300 static constexpr ThreadExitCodeEnum ITK_PROCESS_ABORTED_EXCEPTION =
301 ThreadExitCodeEnum::ITK_PROCESS_ABORTED_EXCEPTION;
302 static constexpr ThreadExitCodeEnum STD_EXCEPTION = ThreadExitCodeEnum::STD_EXCEPTION;
303 static constexpr ThreadExitCodeEnum UNKNOWN = ThreadExitCodeEnum::UNKNOWN;
304#endif
305 };
312 virtual void
314
319 virtual void
321
323 void
325
326#ifndef ITK_FUTURE_LEGACY_REMOVE
327 // `TemplatedThreadingFunctorType` was previously used to declare the `funcP` parameter of `ParallelizeImageRegion`
328 // and `ParallelizeImageRegionRestrictDirection` template member functions.
329 template <unsigned int VDimension>
330 using TemplatedThreadingFunctorType [[deprecated]] = std::function<void(const ImageRegion<VDimension> &)>;
331#endif
332
333 using ThreadingFunctorType = std::function<void(const IndexValueType index[], const SizeValueType size[])>;
334 using ArrayThreadingFunctorType = std::function<void(SizeValueType)>;
335
341 virtual void
343 SizeValueType lastIndexPlus1,
345 ProcessObject * filter);
346
352 template <unsigned int VDimension, typename TFunction>
353 ITK_TEMPLATE_EXPORT void
354 ParallelizeImageRegion(const ImageRegion<VDimension> & requestedRegion, TFunction funcP, ProcessObject * filter)
355 {
356 this->ParallelizeImageRegion(
357 VDimension,
358 requestedRegion.GetIndex().m_InternalArray,
359 requestedRegion.GetSize().m_InternalArray,
360 [&funcP](const IndexValueType index[], const SizeValueType size[]) {
361 ImageRegion<VDimension> region;
362 for (unsigned int d = 0; d < VDimension; ++d)
363 {
364 region.SetIndex(d, index[d]);
365 region.SetSize(d, size[d]);
366 }
367 funcP(region);
368 },
369 filter);
370 }
376 template <unsigned int VDimension, typename TFunction>
377 ITK_TEMPLATE_EXPORT void
378 ParallelizeImageRegionRestrictDirection(unsigned int restrictedDirection,
379 const ImageRegion<VDimension> & requestedRegion,
380 TFunction funcP,
381 ProcessObject * filter)
382 {
383 if constexpr (VDimension <= 1) // Cannot split, no parallelization
384 {
385
386 ProgressReporter progress(filter, 0, requestedRegion.GetNumberOfPixels());
387 funcP(requestedRegion);
388 }
389 else // Can split, parallelize!
390 {
391 constexpr unsigned int SplitDimension = VDimension - 1;
392
393 Index<SplitDimension> splitIndex{};
394 Size<SplitDimension> splitSize{};
395 for (unsigned int splitDimension = 0, dimension = 0; dimension < VDimension; ++dimension)
396 {
397 if (dimension != restrictedDirection)
398 {
399 splitIndex[splitDimension] = requestedRegion.GetIndex(dimension);
400 splitSize[splitDimension] = requestedRegion.GetSize(dimension);
401 ++splitDimension;
402 }
403 }
404
405 this->ParallelizeImageRegion(
406 SplitDimension,
407 splitIndex.m_InternalArray,
408 splitSize.m_InternalArray,
409 [restrictedDirection, &requestedRegion, &funcP](const IndexValueType index[], const SizeValueType size[]) {
410 ImageRegion<VDimension> restrictedRequestedRegion;
411 restrictedRequestedRegion.SetIndex(restrictedDirection, requestedRegion.GetIndex(restrictedDirection));
412 restrictedRequestedRegion.SetSize(restrictedDirection, requestedRegion.GetSize(restrictedDirection));
413 for (unsigned int splitDimension = 0, dimension = 0; dimension < VDimension; ++dimension)
414 {
415 if (dimension != restrictedDirection)
416 {
417 restrictedRequestedRegion.SetIndex(dimension, index[splitDimension]);
418 restrictedRequestedRegion.SetSize(dimension, size[splitDimension]);
419 ++splitDimension;
420 }
421 }
422 funcP(restrictedRequestedRegion);
423 },
424 filter);
425 }
426 }
427
430 virtual void
431 ParallelizeImageRegion(unsigned int dimension,
432 const IndexValueType index[],
433 const SizeValueType size[],
435 ProcessObject * filter);
436
437protected:
440 void
441 PrintSelf(std::ostream & os, Indent indent) const override;
442
444 {
449 };
450
453
455 {
457 unsigned int dimension;
461 };
462
465
467 ThreadIdType m_NumberOfWorkUnits{};
468
477 ThreadIdType m_MaximumNumberOfThreads{};
478
486 SingleMethodProxy(void * arg);
487
489 ThreadFunctionType m_SingleMethod{ nullptr };
490
492 void * m_SingleData{ nullptr };
493
494private:
495 static void
497 static ThreaderEnum
499
501 itkGetGlobalDeclarationMacro(MultiThreaderBaseGlobals, PimplGlobals);
502
503 std::atomic<bool> m_UpdateProgress{ true };
504
505 static MultiThreaderBaseGlobals * m_PimplGlobals;
509 friend class ProcessObject;
510
512 static ThreadIdType
514};
515
516
517} // namespace itk
518
519#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)
itkLegacyMacro(static void SetGlobalDefaultUseThreadPool(const bool GlobalDefaultUseThreadPool);) itkLegacyMacro(static bool GetGlobalDefaultUseThreadPool()
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()
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....
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
unsigned int ThreadIdType
Definition: itkIntTypes.h:102
long IndexValueType
Definition: itkIntTypes.h:93
void(*)(void *) ThreadFunctionType
unsigned long SizeValueType
Definition: itkIntTypes.h:86
itk::ITK_THREAD_RETURN_TYPE ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:69
IndexValueType m_InternalArray[VDimension]
Definition: itkIndex.h:291
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:70
SizeValueType m_InternalArray[VDimension]
Definition: itkSize.h:243