RTK  2.6.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  * 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 
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_TEMPLATE_EXPORT AmsterdamShroudImageFilter
82  : public itk::ImageToImageFilter<TInputImage, itk::Image<double, TInputImage::ImageDimension - 1>>
83 {
84 public:
85  ITK_DISALLOW_COPY_AND_MOVE(AmsterdamShroudImageFilter);
86 
89  using Superclass = itk::ImageToImageFilter<TInputImage, itk::Image<double, TInputImage::ImageDimension - 1>>;
92 
94  using TOutputImage = itk::Image<double, TInputImage::ImageDimension - 1>;
98 
100  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
101  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
102  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
103 
105  itkNewMacro(Self);
106 
112  itkGetMacro(UnsharpMaskSize, unsigned int);
113  itkSetMacro(UnsharpMaskSize, unsigned int);
115 
117  itkGetModifiableObjectMacro(Geometry, GeometryType);
118  itkSetObjectMacro(Geometry, GeometryType);
120 
124  itkGetMacro(Corner1, PointType);
125  itkSetMacro(Corner1, PointType);
126  itkGetMacro(Corner2, PointType);
127  itkSetMacro(Corner2, PointType);
129 
131 #ifdef itkOverrideGetNameOfClassMacro
132  itkOverrideGetNameOfClassMacro(AmsterdamShroudImageFilter);
133 #else
135 #endif
136 
137 
138 protected:
140  ~AmsterdamShroudImageFilter() override = default;
141 
142  void
143  GenerateOutputInformation() override;
144  void
145  GenerateInputRequestedRegion() override;
146  void
147  UpdateUnsharpMaskKernel();
148 
151  void
152  GenerateData() override;
153 
156  virtual void
157  CropOutsideProjectedBox();
158 
159 private:
167 
175  unsigned int m_UnsharpMaskSize{ 17 };
176  GeometryPointer m_Geometry{ nullptr };
177  PointType m_Corner1{ 0. };
178  PointType m_Corner2{ 0. };
179 }; // end of class
180 
181 } // end namespace rtk
182 
183 #ifndef ITK_MANUAL_INSTANTIATION
184 # include "rtkAmsterdamShroudImageFilter.hxx"
185 #endif
186 
187 #endif
Projection geometry for a source and a 2-D flat panel.
#define itkSetMacro(name, type)
Compute the Amsterdam shroud image for respiratory signal extraction.
typename GeometryType::Pointer GeometryPointer