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
84// Define how to print enumeration
85extern ITKCommon_EXPORT std::ostream &
86 operator<<(std::ostream & out, const MultiThreaderBaseEnums::Threader value);
87extern ITKCommon_EXPORT std::ostream &
88 operator<<(std::ostream & out, const MultiThreaderBaseEnums::ThreadExitCode value);
89
103
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
129 virtual void
131 itkGetConstMacro(MaximumNumberOfThreads, ThreadIdType);
137 virtual void
139 itkGetConstMacro(NumberOfWorkUnits, ThreadIdType);
141 virtual void
142 SetUpdateProgress(bool updates);
143 itkGetConstMacro(UpdateProgress, bool);
144
151 static void
153 static ThreadIdType
162 itkLegacyMacro(static void SetGlobalDefaultUseThreadPool(const bool GlobalDefaultUseThreadPool);)
163 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 Single = ThreaderEnum::Single;
173 static constexpr ThreaderEnum Last = ThreaderEnum::Last;
174 static constexpr ThreaderEnum Unknown = ThreaderEnum::Unknown;
175#endif
176
178 static ThreaderEnum
179 ThreaderTypeFromString(std::string threaderString);
180
182 static std::string
184 {
185 switch (threader)
186 {
187 case ThreaderEnum::Platform:
188 return "Platform";
189 case ThreaderEnum::Pool:
190 return "Pool";
191 case ThreaderEnum::TBB:
192 return "TBB";
193 case ThreaderEnum::Single:
194 return "Single";
195 case ThreaderEnum::Unknown:
196 default:
197 return "Unknown";
198 }
199 }
200
213 static void
215 static ThreaderEnum
223 static void
225 static ThreadIdType
228#if !defined(ITK_LEGACY_REMOVE)
232 itkLegacyMacro(virtual void SetNumberOfThreads(ThreadIdType numberOfThreads))
233 {
234 this->SetMaximumNumberOfThreads(numberOfThreads);
235 this->SetNumberOfWorkUnits(this->GetMaximumNumberOfThreads()); // Might be clamped
236 }
237 itkLegacyMacro(virtual ThreadIdType GetNumberOfThreads())
238 {
239 return this->GetNumberOfWorkUnits();
240 }
251 // clang-format off
252ITK_GCC_PRAGMA_DIAG_PUSH()
253ITK_GCC_PRAGMA_DIAG(ignored "-Wattributes")
254INTEL_PRAGMA_WARN_PUSH
255INTEL_SUPPRESS_warning_1292
256 // clang-format on
257# ifdef ITK_LEGACY_SILENT
258 struct ThreadInfoStruct
259# else
260 struct [[deprecated("Use WorkUnitInfo, ThreadInfoStruct is deprecated since ITK 5.0")]] ThreadInfoStruct
261# endif
262 // clang-format off
263INTEL_PRAGMA_WARN_POP
264 // clang-format on
265 {
266 ThreadIdType ThreadID;
267 ThreadIdType NumberOfThreads;
268 void * UserData;
269 ThreadFunctionType ThreadFunction;
270 enum
271 {
272 SUCCESS,
273 ITK_EXCEPTION,
274 ITK_PROCESS_ABORTED_EXCEPTION,
275 STD_EXCEPTION,
276 UNKNOWN
277 } ThreadExitCode;
278 };
279 // clang-format off
280ITK_GCC_PRAGMA_DIAG_POP()
281 // clang-format on
283#endif // ITK_LEGACY_REMOVE
284
292 {
295 void * UserData;
299#if !defined(ITK_LEGACY_REMOVE)
300 using ThreadExitCodeType = ThreadExitCodeEnum;
301 static constexpr ThreadExitCodeEnum SUCCESS = ThreadExitCodeEnum::SUCCESS;
302 static constexpr ThreadExitCodeEnum ITK_EXCEPTION = ThreadExitCodeEnum::ITK_EXCEPTION;
303 static constexpr ThreadExitCodeEnum ITK_PROCESS_ABORTED_EXCEPTION =
304 ThreadExitCodeEnum::ITK_PROCESS_ABORTED_EXCEPTION;
305 static constexpr ThreadExitCodeEnum STD_EXCEPTION = ThreadExitCodeEnum::STD_EXCEPTION;
306 static constexpr ThreadExitCodeEnum UNKNOWN = ThreadExitCodeEnum::UNKNOWN;
307#endif
308 };
309
314 virtual void
316
321 virtual void
323
325 void
327
328#ifndef ITK_FUTURE_LEGACY_REMOVE
329 // `TemplatedThreadingFunctorType` was previously used to declare the `funcP` parameter of `ParallelizeImageRegion`
330 // and `ParallelizeImageRegionRestrictDirection` template member functions.
331 template <unsigned int VDimension>
332 using TemplatedThreadingFunctorType [[deprecated]] = std::function<void(const ImageRegion<VDimension> &)>;
333#endif
334
335 using ThreadingFunctorType = std::function<void(const IndexValueType index[], const SizeValueType size[])>;
336 using ArrayThreadingFunctorType = std::function<void(SizeValueType)>;
337
343 virtual void
345 SizeValueType lastIndexPlus1,
347 ProcessObject * filter);
348
354 template <unsigned int VDimension, typename TFunction>
355 ITK_TEMPLATE_EXPORT void
356 ParallelizeImageRegion(const ImageRegion<VDimension> & requestedRegion, TFunction funcP, ProcessObject * filter)
357 {
359 VDimension,
360 requestedRegion.GetIndex().m_InternalArray,
361 requestedRegion.GetSize().m_InternalArray,
362 [&funcP](const IndexValueType index[], const SizeValueType size[]) {
363 ImageRegion<VDimension> region;
364 for (unsigned int d = 0; d < VDimension; ++d)
365 {
366 region.SetIndex(d, index[d]);
367 region.SetSize(d, size[d]);
368 }
369 funcP(region);
370 },
371 filter);
372 }
373
377 template <unsigned int VDimension, typename TFunction>
378 ITK_TEMPLATE_EXPORT void
379 ParallelizeImageRegionRestrictDirection(unsigned int restrictedDirection,
380 const ImageRegion<VDimension> & requestedRegion,
381 TFunction funcP,
382 ProcessObject * filter)
383 {
384 if constexpr (VDimension <= 1) // Cannot split, no parallelization
385 {
386
387 const ProgressReporter progress(filter, 0, requestedRegion.GetNumberOfPixels());
388 funcP(requestedRegion);
389 }
390 else // Can split, parallelize!
391 {
392 constexpr unsigned int SplitDimension = VDimension - 1;
393
394 Index<SplitDimension> splitIndex{};
395 Size<SplitDimension> splitSize{};
396 for (unsigned int splitDimension = 0, dimension = 0; dimension < VDimension; ++dimension)
397 {
398 if (dimension != restrictedDirection)
399 {
400 splitIndex[splitDimension] = requestedRegion.GetIndex(dimension);
401 splitSize[splitDimension] = requestedRegion.GetSize(dimension);
402 ++splitDimension;
403 }
404 }
405
407 SplitDimension,
408 splitIndex.m_InternalArray,
409 splitSize.m_InternalArray,
410 [restrictedDirection, &requestedRegion, &funcP](const IndexValueType index[], const SizeValueType size[]) {
411 ImageRegion<VDimension> restrictedRequestedRegion;
412 restrictedRequestedRegion.SetIndex(restrictedDirection, requestedRegion.GetIndex(restrictedDirection));
413 restrictedRequestedRegion.SetSize(restrictedDirection, requestedRegion.GetSize(restrictedDirection));
414 for (unsigned int splitDimension = 0, dimension = 0; dimension < VDimension; ++dimension)
415 {
416 if (dimension != restrictedDirection)
417 {
418 restrictedRequestedRegion.SetIndex(dimension, index[splitDimension]);
419 restrictedRequestedRegion.SetSize(dimension, size[splitDimension]);
420 ++splitDimension;
421 }
422 }
423 funcP(restrictedRequestedRegion);
424 },
425 filter);
426 }
427 }
428
431 virtual void
432 ParallelizeImageRegion(unsigned int dimension,
433 const IndexValueType index[],
434 const SizeValueType size[],
436 ProcessObject * filter);
437
438protected:
441 void
442 PrintSelf(std::ostream & os, Indent indent) const override;
443
451
454
463
466
469
479
487 SingleMethodProxy(void * arg);
488
491
493 void * m_SingleData{ nullptr };
494
495private:
496 static void
498 static ThreaderEnum
500
502 itkGetGlobalDeclarationMacro(MultiThreaderBaseGlobals, PimplGlobals);
503
504 std::atomic<bool> m_UpdateProgress{ true };
505
506 static MultiThreaderBaseGlobals * m_PimplGlobals;
510 friend class ProcessObject;
511
513 static ThreadIdType
515};
516
517
518} // namespace itk
519
520#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
void(*)(void *) ThreadFunctionType
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, AnatomicalOrientation::CoordinateEnum value)
ITK_THREAD_RETURN_TYPE ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION
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