RTK  1.4.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  * 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 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>
40  public itk::InPlaceImageFilter<TInputImage,TOutputImage>
41 {
42 public:
48 
51  typedef std::string StringType;
54 
56  itkNewMacro(Self);
57 
60 
62  itkGetConstObjectMacro(GeometricPhantom, GeometricPhantom);
63  itkSetConstObjectMacro(GeometricPhantom, GeometricPhantom);
65 
67  itkSetMacro(ConfigFile, StringType);
68  itkGetMacro(ConfigFile, StringType);
70 
72  itkSetMacro(PhantomScale, VectorType);
73  itkGetMacro(PhantomScale, VectorType);
75 
78  itkSetMacro(OriginOffset, VectorType);
79  itkGetMacro(OriginOffset, VectorType);
81 
84  itkSetMacro(IsForbildConfigFile, bool);
85  itkGetConstMacro(IsForbildConfigFile, bool);
86  itkBooleanMacro(IsForbildConfigFile);
88 
90  itkSetMacro(RotationMatrix, RotationMatrixType);
91  itkGetMacro(RotationMatrix, RotationMatrixType);
93 
94 protected:
96  virtual ~DrawGeometricPhantomImageFilter() ITK_OVERRIDE {}
97 
98  void GenerateData() ITK_OVERRIDE;
99 
100 private:
101  DrawGeometricPhantomImageFilter(const Self&); //purposely not implemented
102  void operator=(const Self&); //purposely not implemented
103 
104  GeometricPhantomConstPointer m_GeometricPhantom;
105  StringType m_ConfigFile;
106  VectorType m_PhantomScale;
107  VectorType m_OriginOffset;
109  RotationMatrixType m_RotationMatrix;
110 };
111 
112 } // end namespace rtk
113 
114 #ifndef ITK_MANUAL_INSTANTIATION
115 #include "rtkDrawGeometricPhantomImageFilter.hxx"
116 #endif
117 
118 #endif
Container for a geometric phantom, i.e., a set of ConvexShapes.
GeometricPhantom::ConstPointer GeometricPhantomConstPointer
Draws a GeometricPhantom in a 3D image.
itk::InPlaceImageFilter< TInputImage, TOutputImage > Superclass
#define itkSetMacro(name, type)