ITK  5.2.0
Insight Toolkit
itkThreadedImageRegionPartitioner.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  * http://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 #ifndef itkThreadedImageRegionPartitioner_h
19 #define itkThreadedImageRegionPartitioner_h
20 
22 #include "itkImageRegion.h"
24 
25 namespace itk
26 {
27 
44 template <unsigned int VDimension>
45 class ITK_TEMPLATE_EXPORT ThreadedImageRegionPartitioner : public ThreadedDomainPartitioner<ImageRegion<VDimension>>
46 {
47 public:
48  ITK_DISALLOW_COPY_AND_MOVE(ThreadedImageRegionPartitioner);
49 
55 
57  itkNewMacro(Self);
58 
61 
63  using DomainType = typename Superclass::DomainType;
64 
66  static constexpr unsigned int ImageDimension = VDimension;
67  using ImageRegionType = typename Self::DomainType;
70 
84  PartitionDomain(const ThreadIdType threadId,
85  const ThreadIdType requestedTotal,
86  const DomainType & completeRegion,
87  DomainType & subRegion) const override;
88 
89 protected:
91  ~ThreadedImageRegionPartitioner() override = default;
92 
94 
95 private:
97 };
98 
99 } // end namespace itk
100 
101 #ifndef ITK_MANUAL_INSTANTIATION
102 # include "itkThreadedImageRegionPartitioner.hxx"
103 #endif
104 
105 #endif
itk::ImageRegionSplitterSlowDimension
Divide an image region along the slowest dimension.
Definition: itkImageRegionSplitterSlowDimension.h:46
itk::ThreadedImageRegionPartitioner::SizeType
typename Self::DomainType::SizeType SizeType
Definition: itkThreadedImageRegionPartitioner.h:68
itk::ThreadedImageRegionPartitioner
Class for partitioning of an ImageRegion.
Definition: itkThreadedImageRegionPartitioner.h:45
itk::ThreadedImageRegionPartitioner::m_ImageRegionSplitter
ImageRegionSplitterType::Pointer m_ImageRegionSplitter
Definition: itkThreadedImageRegionPartitioner.h:96
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itkThreadedDomainPartitioner.h
itkImageRegion.h
itk::ThreadIdType
unsigned int ThreadIdType
Definition: itkIntTypes.h:99
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itkImageRegionSplitterSlowDimension.h
itk::ThreadedImageRegionPartitioner::DomainType
typename Superclass::DomainType DomainType
Definition: itkThreadedImageRegionPartitioner.h:63
itk::ThreadedImageRegionPartitioner::IndexType
typename Self::DomainType::IndexType IndexType
Definition: itkThreadedImageRegionPartitioner.h:69
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ThreadedDomainPartitioner
Virtual base class for partitioning a domain into subsets to be processed per thread when parallel pr...
Definition: itkThreadedDomainPartitioner.h:47
itk::ThreadedImageRegionPartitioner::ImageRegionType
typename Self::DomainType ImageRegionType
Definition: itkThreadedImageRegionPartitioner.h:67