RTK  2.6.0
Reconstruction Toolkit
rtkDrawGeometricPhantomImageFilter.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 rtkDrawGeometricPhantomImageFilter_h
20 #define rtkDrawGeometricPhantomImageFilter_h
21 
22 #include <itkInPlaceImageFilter.h>
23 #include <itkAddImageFilter.h>
24 #include "rtkGeometricPhantom.h"
25 
26 namespace rtk
27 {
28 
38 template <class TInputImage, class TOutputImage>
39 class ITK_TEMPLATE_EXPORT DrawGeometricPhantomImageFilter : public itk::InPlaceImageFilter<TInputImage, TOutputImage>
40 {
41 public:
42  ITK_DISALLOW_COPY_AND_MOVE(DrawGeometricPhantomImageFilter);
43 
49 
52  using StringType = std::string;
56 
58  itkNewMacro(Self);
59 
61  itkOverrideGetNameOfClassMacro(DrawGeometricPhantomImageFilter);
62 
64  itkGetConstObjectMacro(GeometricPhantom, GeometricPhantom);
65  itkSetConstObjectMacro(GeometricPhantom, GeometricPhantom);
67 
69  itkSetMacro(ConfigFile, StringType);
70  itkGetMacro(ConfigFile, StringType);
72 
74  itkSetMacro(PhantomScale, VectorType);
75  itkGetMacro(PhantomScale, VectorType);
77 
79  virtual void
81  {
82  SetPhantomScale(VectorType(_arg));
83  }
84 
87  itkSetMacro(OriginOffset, VectorType);
88  itkGetMacro(OriginOffset, VectorType);
90 
93  itkSetMacro(IsForbildConfigFile, bool);
94  itkGetConstMacro(IsForbildConfigFile, bool);
95  itkBooleanMacro(IsForbildConfigFile);
97 
99  itkSetMacro(RotationMatrix, RotationMatrixType);
100  itkGetMacro(RotationMatrix, RotationMatrixType);
102 
105  void
106  AddClipPlane(const VectorType & dir, const ScalarType & pos);
107  void
108  SetClipPlanes(const std::vector<VectorType> & dir, const std::vector<ScalarType> & pos);
110 
111 protected:
113  ~DrawGeometricPhantomImageFilter() override = default;
114 
115  void
116  GenerateData() override;
117 
118 private:
121  VectorType m_PhantomScale{ 1. };
122  VectorType m_OriginOffset{ 0. };
123  bool m_IsForbildConfigFile{ false };
125  std::vector<VectorType> m_PlaneDirections;
126  std::vector<ScalarType> m_PlanePositions;
127 };
128 
129 } // end namespace rtk
130 
131 #ifndef ITK_MANUAL_INSTANTIATION
132 # include "rtkDrawGeometricPhantomImageFilter.hxx"
133 #endif
134 
135 #endif
Container for a geometric phantom, i.e., a set of ConvexShapes.
itk::Vector< ScalarType, Dimension > VectorType
virtual void SetPhantomScale(const ScalarType _arg)
itk::SmartPointer< const Self > ConstPointer
Draws a GeometricPhantom in a 3D image.
#define itkSetMacro(name, type)
itk::Matrix< ScalarType, Dimension, Dimension > RotationMatrixType