RTK  1.4.0
Reconstruction Toolkit
rtkAverageOutOfROIImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright RTK 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 rtkAverageOutOfROIImageFilter_h
19 #define rtkAverageOutOfROIImageFilter_h
20 
21 #include "itkInPlaceImageFilter.h"
22 
24 
25 #include "rtkMacro.h"
26 
27 namespace rtk
28 {
50 template< class TInputImage,
51  class TROI = itk::Image< typename TInputImage::PixelType, TInputImage::ImageDimension -1 > >
52 
53 class AverageOutOfROIImageFilter : public itk::InPlaceImageFilter<TInputImage, TInputImage>
54 {
55 public:
60  typedef itk::Image<typename TInputImage::PixelType, TInputImage::ImageDimension - 1> LowerDimImage;
61 
63  itkNewMacro(Self)
64 
65 
67 
69  void SetROI(const TROI* Map);
70 
71 protected:
72  AverageOutOfROIImageFilter();
73  virtual ~AverageOutOfROIImageFilter() ITK_OVERRIDE {}
74 
75  typename TROI::Pointer GetROI();
76 
77  void GenerateOutputInformation() ITK_OVERRIDE;
78  void GenerateInputRequestedRegion() ITK_OVERRIDE;
79 
81  void ThreadedGenerateData(const typename TInputImage::RegionType& outputRegionForThread, itk::ThreadIdType itkNotUsed(threadId)) ITK_OVERRIDE;
82 
84  const itk::ImageRegionSplitterBase* GetImageRegionSplitter(void) const ITK_OVERRIDE;
86 
87 private:
88  AverageOutOfROIImageFilter(const Self &); //purposely not implemented
89  void operator=(const Self &); //purposely not implemented
90 
91 };
92 } //namespace RTK
93 
94 
95 #ifndef ITK_MANUAL_INSTANTIATION
96 #include "rtkAverageOutOfROIImageFilter.hxx"
97 #endif
98 
99 #endif
itk::ImageRegionSplitterDirection::Pointer m_Splitter
const itk::ImageRegionSplitterBase * GetImageRegionSplitter(void) const override
void GenerateOutputInformation() override
itk::Image< typename TInputImage::PixelType, TInputImage::ImageDimension-1 > LowerDimImage
void GenerateInputRequestedRegion() override
void ThreadedGenerateData(const typename TInputImage::RegionType &outputRegionForThread, itk::ThreadIdType) override
void operator=(const Self &)
unsigned int ThreadIdType
itk::ImageToImageFilter< TInputImage, TInputImage > Superclass
void SetROI(const TROI *Map)
Averages along the last dimension if the pixel is outside ROI.