ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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
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
102
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
128 virtual void
130 itkGetConstMacro(MaximumNumberOfThreads, ThreadIdType);
136 virtual void
138 itkGetConstMacro(NumberOfWorkUnits, ThreadIdType);
140 virtual void
141 SetUpdateProgress(bool updates);
142 itkGetConstMacro(UpdateProgress, bool);
143
150 static void
152 static ThreadIdType
161 itkLegacyMacro(static void SetGlobalDefaultUseThreadPool(const bool GlobalDefaultUseThreadPool);)
162 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 }
196
209 static void
211 static ThreaderEnum
219 static void
221 static ThreadIdType
224#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())
234 {
235 return this->GetNumberOfWorkUnits();
236 }
247 // clang-format off
248ITK_GCC_PRAGMA_DIAG_PUSH()
249ITK_GCC_PRAGMA_DIAG(ignored "-Wattributes")
250INTEL_PRAGMA_WARN_PUSH
251INTEL_SUPPRESS_warning_1292
252 // clang-format on
253# ifdef ITK_LEGACY_SILENT
254 struct ThreadInfoStruct
255# else
256 struct [[deprecated("Use WorkUnitInfo, ThreadInfoStruct is deprecated since ITK 5.0")]] ThreadInfoStruct
257# endif
258 // clang-format off
259INTEL_PRAGMA_WARN_POP
260 // clang-format on
261 {
262 ThreadIdType ThreadID;
263 ThreadIdType NumberOfThreads;
264 void * UserData;
265 ThreadFunctionType ThreadFunction;
266 enum
267 {
268 SUCCESS,
269 ITK_EXCEPTION,
270 ITK_PROCESS_ABORTED_EXCEPTION,
271 STD_EXCEPTION,
272 UNKNOWN
273 } ThreadExitCode;
274 };
275 // clang-format off
276ITK_GCC_PRAGMA_DIAG_POP()
277 // clang-format on
279#endif // ITK_LEGACY_REMOVE
280
288 {
291 void * UserData;
295#if !defined(ITK_LEGACY_REMOVE)
296 using ThreadExitCodeType = ThreadExitCodeEnum;
297 static constexpr ThreadExitCodeEnum SUCCESS = ThreadExitCodeEnum::SUCCESS;
298 static constexpr ThreadExitCodeEnum ITK_EXCEPTION = ThreadExitCodeEnum::ITK_EXCEPTION;
299 static constexpr ThreadExitCodeEnum ITK_PROCESS_ABORTED_EXCEPTION =
300 ThreadExitCodeEnum::ITK_PROCESS_ABORTED_EXCEPTION;
301 static constexpr ThreadExitCodeEnum STD_EXCEPTION = ThreadExitCodeEnum::STD_EXCEPTION;
302 static constexpr ThreadExitCodeEnum UNKNOWN = ThreadExitCodeEnum::UNKNOWN;
303#endif
304 };
305
310 virtual void
312
317 virtual void
319
321 void
323
324#ifndef ITK_FUTURE_LEGACY_REMOVE
325 // `TemplatedThreadingFunctorType` was previously used to declare the `funcP` parameter of `ParallelizeImageRegion`
326 // and `ParallelizeImageRegionRestrictDirection` template member functions.
327 template <unsigned int VDimension>
328 using TemplatedThreadingFunctorType [[deprecated]] = std::function<void(const ImageRegion<VDimension> &)>;
329#endif
330
331 using ThreadingFunctorType = std::function<void(const IndexValueType index[], const SizeValueType size[])>;
332 using ArrayThreadingFunctorType = std::function<void(SizeValueType)>;
333
339 virtual void
341 SizeValueType lastIndexPlus1,
343 ProcessObject * filter);
344
350 template <unsigned int VDimension, typename TFunction>
351 ITK_TEMPLATE_EXPORT void
352 ParallelizeImageRegion(const ImageRegion<VDimension> & requestedRegion, TFunction funcP, ProcessObject * filter)
353 {
355 VDimension,
356 requestedRegion.GetIndex().m_InternalArray,
357 requestedRegion.GetSize().m_InternalArray,
358 [&funcP](const IndexValueType index[], const SizeValueType size[]) {
359 ImageRegion<VDimension> region;
360 for (unsigned int d = 0; d < VDimension; ++d)
361 {
362 region.SetIndex(d, index[d]);
363 region.SetSize(d, size[d]);
364 }
365 funcP(region);
366 },
367 filter);
368 }
369
373 template <unsigned int VDimension, typename TFunction>
374 ITK_TEMPLATE_EXPORT void
375 ParallelizeImageRegionRestrictDirection(unsigned int restrictedDirection,
376 const ImageRegion<VDimension> & requestedRegion,
377 TFunction funcP,
378 ProcessObject * filter)
379 {
380 if constexpr (VDimension <= 1) // Cannot split, no parallelization
381 {
382
383 const ProgressReporter progress(filter, 0, requestedRegion.GetNumberOfPixels());
384 funcP(requestedRegion);
385 }
386 else // Can split, parallelize!
387 {
388 constexpr unsigned int SplitDimension = VDimension - 1;
389
390 Index<SplitDimension> splitIndex{};
391 Size<SplitDimension> splitSize{};
392 for (unsigned int splitDimension = 0, dimension = 0; dimension < VDimension; ++dimension)
393 {
394 if (dimension != restrictedDirection)
395 {
396 splitIndex[splitDimension] = requestedRegion.GetIndex(dimension);
397 splitSize[splitDimension] = requestedRegion.GetSize(dimension);
398 ++splitDimension;
399 }
400 }
401
403 SplitDimension,
404 splitIndex.m_InternalArray,
405 splitSize.m_InternalArray,
406 [restrictedDirection, &requestedRegion, &funcP](const IndexValueType index[], const SizeValueType size[]) {
407 ImageRegion<VDimension> restrictedRequestedRegion;
408 restrictedRequestedRegion.SetIndex(restrictedDirection, requestedRegion.GetIndex(restrictedDirection));
409 restrictedRequestedRegion.SetSize(restrictedDirection, requestedRegion.GetSize(restrictedDirection));
410 for (unsigned int splitDimension = 0, dimension = 0; dimension < VDimension; ++dimension)
411 {
412 if (dimension != restrictedDirection)
413 {
414 restrictedRequestedRegion.SetIndex(dimension, index[splitDimension]);
415 restrictedRequestedRegion.SetSize(dimension, size[splitDimension]);
416 ++splitDimension;
417 }
418 }
419 funcP(restrictedRequestedRegion);
420 },
421 filter);
422 }
423 }
424
427 virtual void
428 ParallelizeImageRegion(unsigned int dimension,
429 const IndexValueType index[],
430 const SizeValueType size[],
432 ProcessObject * filter);
433
434protected:
437 void
438 PrintSelf(std::ostream & os, Indent indent) const override;
439
447
450
459
462
465
475
483 SingleMethodProxy(void * arg);
484
487
489 void * m_SingleData{ nullptr };
490
491private:
492 static void
494 static ThreaderEnum
496
498 itkGetGlobalDeclarationMacro(MultiThreaderBaseGlobals, PimplGlobals);
499
500 std::atomic<bool> m_UpdateProgress{ true };
501
502 static MultiThreaderBaseGlobals * m_PimplGlobals;
506 friend class ProcessObject;
507
509 static ThreadIdType
511};
512
513
514} // namespace itk
515
516#endif
An image region represents a structured region of data.
SizeValueType GetNumberOfPixels() const
const IndexType & GetIndex() const
const SizeType & GetSize() const
Control indentation during Print() invocation.
Definition itkIndent.h:50
enums for MultiThreaderBase
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ParallelizeArrayHelper(void *arg)
static bool GetGlobalDefaultUseThreadPool()
ITK_TEMPLATE_EXPORT void ParallelizeImageRegion(const ImageRegion< VDimension > &requestedRegion, TFunction funcP, ProcessObject *filter)
MultiThreaderBaseEnums::Threader ThreaderEnum
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 void SetGlobalDefaultUseThreadPool(const bool GlobalDefaultUseThreadPool)
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)
SmartPointer< Self > Pointer
static MultiThreaderBaseGlobals * m_PimplGlobals
static ThreaderEnum ThreaderTypeFromString(std::string threaderString)
ThreadFunctionType m_SingleMethod
static ThreadIdType GetGlobalMaximumNumberOfThreads()
virtual void SetNumberOfWorkUnits(ThreadIdType numberOfWorkUnits)
SmartPointer< const Self > ConstPointer
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 ThreadIdType GetMaximumNumberOfThreads() const
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)
std::atomic< bool > m_UpdateProgress
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Implements progress tracking for a filter.
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned int ThreadIdType
unsigned long SizeValueType
Definition itkIntTypes.h:86
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
itk::ITK_THREAD_RETURN_TYPE ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
void(*)(void *) ThreadFunctionType
long IndexValueType
Definition itkIntTypes.h:93
Represent a n-dimensional index in a n-dimensional image.
Definition itkIndex.h:69
IndexValueType m_InternalArray[VDimension]
Definition itkIndex.h:289
MultiThreaderBaseEnums::ThreadExitCode ThreadExitCodeEnum
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition itkSize.h:70
SizeValueType m_InternalArray[VDimension]
Definition itkSize.h:242