ITK  6.0.0
Insight Toolkit
itkAccumulateImageFilter.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 #ifndef itkAccumulateImageFilter_h
19 #define itkAccumulateImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 
23 namespace itk
24 {
53 template <typename TInputImage, typename TOutputImage>
54 class ITK_TEMPLATE_EXPORT AccumulateImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
55 {
56 public:
57  ITK_DISALLOW_COPY_AND_MOVE(AccumulateImageFilter);
58 
64 
66  itkNewMacro(Self);
67 
69  itkOverrideGetNameOfClassMacro(AccumulateImageFilter);
70 
72  using InputImageType = TInputImage;
75  using InputImagePixelType = typename InputImageType::PixelType;
76 
77  using OutputImageType = TOutputImage;
80  using OutputImagePixelType = typename OutputImageType::PixelType;
81 
83  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
84  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
85 
88 #ifdef ITK_USE_CONCEPT_CHECKING
89  // Begin concept checking
91  // End concept checking
92 #endif
93 
97  itkGetConstMacro(AccumulateDimension, unsigned int);
98  itkSetMacro(AccumulateDimension, unsigned int);
106  itkSetMacro(Average, bool);
107  itkGetConstMacro(Average, bool);
108  itkBooleanMacro(Average);
111 protected:
113  ~AccumulateImageFilter() override = default;
114  void
115  PrintSelf(std::ostream & os, Indent indent) const override;
116 
118  void
119  GenerateOutputInformation() override;
120 
122  void
123  GenerateInputRequestedRegion() override;
124 
129  void
130  GenerateData() override;
131 
132 private:
133  unsigned int m_AccumulateDimension{};
134  bool m_Average{};
135 };
136 } // end namespace itk
137 
138 #ifndef ITK_MANUAL_INSTANTIATION
139 # include "itkAccumulateImageFilter.hxx"
140 #endif
141 
142 #endif
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::Concept::SameDimension
Definition: itkConceptChecking.h:696
itk::ImageToImageFilter::InputImagePixelType
typename InputImageType::PixelType InputImagePixelType
Definition: itkImageToImageFilter.h:133
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::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::AccumulateImageFilter
Implements an accumulation of an image along a selected direction.
Definition: itkAccumulateImageFilter.h:54
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itkImageToImageFilter.h
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::ImageSource::OutputImagePixelType
typename OutputImageType::PixelType OutputImagePixelType
Definition: itkImageSource.h:93
itk::ImageToImageFilter::InputImageRegionType
typename InputImageType::RegionType InputImageRegionType
Definition: itkImageToImageFilter.h:132
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90