RTK  1.4.0
Reconstruction Toolkit
rtkAmsterdamShroudImageFilter.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 
19 #ifndef rtkAmsterdamShroudImageFilter_h
20 #define rtkAmsterdamShroudImageFilter_h
21 
22 #include <itkImageToImageFilter.h>
24 #include <itkMultiplyImageFilter.h>
28 #include <itkSubtractImageFilter.h>
30 #include <itkCropImageFilter.h>
32 
33 namespace rtk
34 {
35 
80 template<class TInputImage>
81 class ITK_EXPORT AmsterdamShroudImageFilter :
82  public itk::ImageToImageFilter<TInputImage, itk::Image<double, TInputImage::ImageDimension-1> >
83 {
84 public:
87  typedef itk::ImageToImageFilter<TInputImage,
88  itk::Image<double, TInputImage::ImageDimension-1> > Superclass;
91 
93  typedef itk::Image<double, TInputImage::ImageDimension-1> TOutputImage;
97 
99  itkStaticConstMacro(InputImageDimension, unsigned int,
100  TInputImage::ImageDimension);
101  itkStaticConstMacro(OutputImageDimension, unsigned int,
102  TOutputImage::ImageDimension);
103  itkStaticConstMacro(ImageDimension, unsigned int,
104  TOutputImage::ImageDimension);
106 
108  itkNewMacro(Self);
109 
115  itkGetMacro(UnsharpMaskSize, unsigned int);
116  itkSetMacro(UnsharpMaskSize, unsigned int);
118 
120  itkGetModifiableObjectMacro(Geometry, GeometryType);
121  itkSetObjectMacro(Geometry, GeometryType);
123 
127  itkGetMacro(Corner1, PointType);
128  itkSetMacro(Corner1, PointType);
129  itkGetMacro(Corner2, PointType);
130  itkSetMacro(Corner2, PointType);
132 
135 protected:
137  virtual ~AmsterdamShroudImageFilter() ITK_OVERRIDE {}
139 
140  void GenerateOutputInformation() ITK_OVERRIDE;
141  void GenerateInputRequestedRegion() ITK_OVERRIDE;
142  void UpdateUnsharpMaskKernel();
143 
146  void GenerateData() ITK_OVERRIDE;
147 
150  virtual void CropOutsideProjectedBox();
151 
152 private:
153  AmsterdamShroudImageFilter(const Self&); //purposely not implemented
154  void operator=(const Self&); //purposely not implemented
155 
163 
171  unsigned int m_UnsharpMaskSize;
172  GeometryPointer m_Geometry;
173  PointType m_Corner1;
174  PointType m_Corner2;
175 }; // end of class
176 
177 } // end namespace rtk
178 
179 #ifndef ITK_MANUAL_INSTANTIATION
180 #include "rtkAmsterdamShroudImageFilter.hxx"
181 #endif
182 
183 #endif
itk::PermuteAxesImageFilter< TOutputImage > PermuteType
itk::ThresholdImageFilter< TInputImage > ThresholdType
itk::RecursiveGaussianImageFilter< TInputImage, TInputImage > DerivativeType
itk::SubtractImageFilter< TOutputImage, TOutputImage > SubtractType
itk::ConvolutionImageFilter< TOutputImage, TOutputImage > ConvolutionType
itk::MultiplyImageFilter< TInputImage, TInputImage, TInputImage > NegativeType
Projection geometry for a source and a 2-D flat panel.
itk::Image< double, TInputImage::ImageDimension-1 > TOutputImage
Compute the Amsterdam shroud image for respiratory signal extraction.
itk::SumProjectionImageFilter< TInputImage, TOutputImage > SumType
itk::SmartPointer< const Self > ConstPointer
rtk::ThreeDCircularProjectionGeometry GeometryType
itk::ImageToImageFilter< TInputImage, itk::Image< double, TInputImage::ImageDimension-1 > > Superclass
#define itkSetMacro(name, type)