RTK  2.6.0
Reconstruction Toolkit
rtkProjectGeometricPhantomImageFilter.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 rtkProjectGeometricPhantomImageFilter_h
20 #define rtkProjectGeometricPhantomImageFilter_h
21 
22 #include <itkInPlaceImageFilter.h>
23 #include <itkAddImageFilter.h>
24 #include "rtkGeometricPhantom.h"
26 
27 namespace rtk
28 {
29 
39 template <class TInputImage, class TOutputImage>
40 class ITK_TEMPLATE_EXPORT ProjectGeometricPhantomImageFilter : public itk::InPlaceImageFilter<TInputImage, TOutputImage>
41 {
42 public:
43  ITK_DISALLOW_COPY_AND_MOVE(ProjectGeometricPhantomImageFilter);
44 
50 
55  using StringType = std::string;
59 
61  itkNewMacro(Self);
62 
64  itkOverrideGetNameOfClassMacro(ProjectGeometricPhantomImageFilter);
65 
67  itkGetConstObjectMacro(GeometricPhantom, GeometricPhantom);
68  itkSetConstObjectMacro(GeometricPhantom, GeometricPhantom);
70 
72  itkGetConstObjectMacro(Geometry, GeometryType);
73  itkSetConstObjectMacro(Geometry, GeometryType);
75 
77  itkSetMacro(ConfigFile, StringType);
78  itkGetMacro(ConfigFile, StringType);
80 
82  itkSetMacro(PhantomScale, VectorType);
83  itkGetMacro(PhantomScale, VectorType);
85 
87  virtual void
89  {
90  SetPhantomScale(VectorType(_arg));
91  }
92 
95  itkSetMacro(OriginOffset, VectorType);
96  itkGetMacro(OriginOffset, VectorType);
98 
101  itkSetMacro(IsForbildConfigFile, bool);
102  itkGetConstMacro(IsForbildConfigFile, bool);
103  itkBooleanMacro(IsForbildConfigFile);
105 
107  itkSetMacro(RotationMatrix, RotationMatrixType);
108  itkGetMacro(RotationMatrix, RotationMatrixType);
110 
113  void
114  AddClipPlane(const VectorType & dir, const ScalarType & pos);
115  void
116  SetClipPlanes(const std::vector<VectorType> & dir, const std::vector<ScalarType> & pos);
118 
119 protected:
121  ~ProjectGeometricPhantomImageFilter() override = default;
122 
124  void
125  VerifyPreconditions() const override;
126 
127  void
128  GenerateData() override;
129 
130 private:
134  VectorType m_PhantomScale{ 1. };
135  VectorType m_OriginOffset{ 0. };
136  bool m_IsForbildConfigFile{ false };
138  std::vector<VectorType> m_PlaneDirections;
139  std::vector<ScalarType> m_PlanePositions;
140 };
141 
142 } // end namespace rtk
143 
144 #ifndef ITK_MANUAL_INSTANTIATION
145 # include "rtkProjectGeometricPhantomImageFilter.hxx"
146 #endif
147 
148 #endif
Container for a geometric phantom, i.e., a set of ConvexShapes.
itk::Vector< ScalarType, Dimension > VectorType
Analytical projection a GeometricPhantom.
itk::SmartPointer< const Self > ConstPointer
Projection geometry for a source and a 2-D flat panel.
#define itkSetMacro(name, type)
itk::Matrix< ScalarType, Dimension, Dimension > RotationMatrixType