RTK  2.0.1
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:
85  ITK_DISALLOW_COPY_AND_ASSIGN(AmsterdamShroudImageFilter);
86 
89  using Superclass = itk::ImageToImageFilter<TInputImage,
90  itk::Image<double, TInputImage::ImageDimension-1> >;
93 
95  using TOutputImage = itk::Image<double, TInputImage::ImageDimension-1>;
99 
101  itkStaticConstMacro(InputImageDimension, unsigned int,
102  TInputImage::ImageDimension);
103  itkStaticConstMacro(OutputImageDimension, unsigned int,
104  TOutputImage::ImageDimension);
105  itkStaticConstMacro(ImageDimension, unsigned int,
106  TOutputImage::ImageDimension);
108 
110  itkNewMacro(Self);
111 
117  itkGetMacro(UnsharpMaskSize, unsigned int);
118  itkSetMacro(UnsharpMaskSize, unsigned int);
120 
122  itkGetModifiableObjectMacro(Geometry, GeometryType);
123  itkSetObjectMacro(Geometry, GeometryType);
125 
129  itkGetMacro(Corner1, PointType);
130  itkSetMacro(Corner1, PointType);
131  itkGetMacro(Corner2, PointType);
132  itkSetMacro(Corner2, PointType);
134 
137 protected:
139  ~AmsterdamShroudImageFilter() override = default;
141 
142  void GenerateOutputInformation() override;
143  void GenerateInputRequestedRegion() override;
144  void UpdateUnsharpMaskKernel();
145 
148  void GenerateData() override;
149 
152  virtual void CropOutsideProjectedBox();
153 
154 private:
162 
170  unsigned int m_UnsharpMaskSize{17};
171  GeometryPointer m_Geometry{nullptr};
172  PointType m_Corner1{0.};
173  PointType m_Corner2{0.};
174 }; // end of class
175 
176 } // end namespace rtk
177 
178 #ifndef ITK_MANUAL_INSTANTIATION
179 #include "rtkAmsterdamShroudImageFilter.hxx"
180 #endif
181 
182 #endif
Projection geometry for a source and a 2-D flat panel.
Compute the Amsterdam shroud image for respiratory signal extraction.
typename GeometryType::Pointer GeometryPointer
#define itkSetMacro(name, type)