ITK  5.2.0
Insight Toolkit
itkRegionOfInterestImageFilter.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 itkRegionOfInterestImageFilter_h
19 #define itkRegionOfInterestImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkSmartPointer.h"
23 
24 namespace itk
25 {
53 template <typename TInputImage, typename TOutputImage>
54 class ITK_TEMPLATE_EXPORT RegionOfInterestImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
55 {
56 public:
57  ITK_DISALLOW_COPY_AND_MOVE(RegionOfInterestImageFilter);
58 
64  using InputImageRegionType = typename Superclass::InputImageRegionType;
65 
67  itkNewMacro(Self);
68 
71 
75  using SizeType = typename TInputImage::SizeType;
76 
78  using OutputImagePixelType = typename TOutputImage::PixelType;
79  using InputImagePixelType = typename TInputImage::PixelType;
80 
82  itkSetMacro(RegionOfInterest, RegionType);
83  itkGetConstMacro(RegionOfInterest, RegionType);
85 
87  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
88  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
89 
90 #ifdef ITK_USE_CONCEPT_CHECKING
91  // Begin concept checking
94  // End concept checking
95 #endif
96 
97 protected:
99  ~RegionOfInterestImageFilter() override = default;
100  void
101  PrintSelf(std::ostream & os, Indent indent) const override;
102 
103  void
104  GenerateInputRequestedRegion() override;
105 
106  void
107  EnlargeOutputRequestedRegion(DataObject * output) override;
108 
116  void
117  GenerateOutputInformation() override;
118 
128  void
129  DynamicThreadedGenerateData(const RegionType & outputRegionForThread) override;
130 
131 private:
133 };
134 } // end namespace itk
135 
136 #ifndef ITK_MANUAL_INSTANTIATION
137 # include "itkRegionOfInterestImageFilter.hxx"
138 #endif
139 
140 #endif
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::Concept::SameDimension
Definition: itkConceptChecking.h:694
itk::ImageToImageFilter::InputImagePixelType
typename InputImageType::PixelType InputImagePixelType
Definition: itkImageToImageFilter.h:133
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::RegionOfInterestImageFilter::IndexType
typename TInputImage::IndexType IndexType
Definition: itkRegionOfInterestImageFilter.h:74
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::RegionOfInterestImageFilter::SizeType
typename TInputImage::SizeType SizeType
Definition: itkRegionOfInterestImageFilter.h:75
itkImageToImageFilter.h
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::Concept::Convertible
Definition: itkConceptChecking.h:216
itkSmartPointer.h
itk::RegionOfInterestImageFilter::RegionType
typename TInputImage::RegionType RegionType
Definition: itkRegionOfInterestImageFilter.h:73
itk::RegionOfInterestImageFilter
Extract a region of interest from the input image.
Definition: itkRegionOfInterestImageFilter.h:54
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::DataObject
Base class for all data objects in ITK.
Definition: itkDataObject.h:293
itk::RegionOfInterestImageFilter::m_RegionOfInterest
RegionType m_RegionOfInterest
Definition: itkRegionOfInterestImageFilter.h:132