ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkProjectionImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 itkProjectionImageFilter_h
19 #define itkProjectionImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
56 template< typename TInputImage, typename TOutputImage, typename TAccumulator >
57 class ITK_TEMPLATE_EXPORT ProjectionImageFilter:
58  public ImageToImageFilter< TInputImage, TOutputImage >
59 {
60 public:
61  ITK_DISALLOW_COPY_AND_ASSIGN(ProjectionImageFilter);
62 
68 
70  itkNewMacro(Self);
71 
74 
76  using InputImageType = TInputImage;
77  using InputImagePointer = typename InputImageType::Pointer;
79  using InputImagePixelType = typename InputImageType::PixelType;
80  using OutputImageType = TOutputImage;
81  using OutputImagePointer = typename OutputImageType::Pointer;
83  using OutputImagePixelType = typename OutputImageType::PixelType;
84 
85  using AccumulatorType = TAccumulator;
86 
88  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
89  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
90 
93 #ifdef ITK_USE_CONCEPT_CHECKING
94  // Begin concept checking
95  itkConceptMacro( ImageDimensionCheck,
97  Self::InputImageDimension,
98  Self::OutputImageDimension > ) );
99  // End concept checking
100 #endif
101 
104  itkSetMacro(ProjectionDimension, unsigned int);
105  itkGetConstReferenceMacro(ProjectionDimension, unsigned int);
107 
108 protected:
110  ~ProjectionImageFilter() override = default;
111  void PrintSelf(std::ostream & os, Indent indent) const override;
112 
114  void GenerateOutputInformation() override;
115 
117  void GenerateInputRequestedRegion() override;
118 
119  void DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
120 
121 
122  virtual AccumulatorType NewAccumulator( SizeValueType ) const;
123 
124 private:
125  unsigned int m_ProjectionDimension;
126 };
127 } // end namespace itk
128 
129 #ifndef ITK_MANUAL_INSTANTIATION
130 #include "itkProjectionImageFilter.hxx"
131 #endif
132 
133 #endif
typename OutputImageType::Pointer OutputImagePointer
unsigned long SizeValueType
Definition: itkIntTypes.h:83
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
Implements an accumulation of an image along a selected direction.
typename OutputImageType::PixelType OutputImagePixelType
typename InputImageType::PixelType InputImagePixelType
typename InputImageType::Pointer InputImagePointer
typename OutputImageType::RegionType OutputImageRegionType
TOutputImage OutputImageType
typename InputImageType::RegionType InputImageRegionType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
#define itkConceptMacro(name, concept)