RTK  2.0.1
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:
43  ITK_DISALLOW_COPY_AND_ASSIGN(DrawGeometricPhantomImageFilter);
44 
50 
53  using StringType = std::string;
56 
58  itkNewMacro(Self);
59 
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 
80  itkSetMacro(OriginOffset, VectorType);
81  itkGetMacro(OriginOffset, VectorType);
83 
86  itkSetMacro(IsForbildConfigFile, bool);
87  itkGetConstMacro(IsForbildConfigFile, bool);
88  itkBooleanMacro(IsForbildConfigFile);
90 
92  itkSetMacro(RotationMatrix, RotationMatrixType);
93  itkGetMacro(RotationMatrix, RotationMatrixType);
95 
96 protected:
98  ~DrawGeometricPhantomImageFilter() override = default;
99 
100  void GenerateData() override;
101 
102 private:
109 };
110 
111 } // end namespace rtk
112 
113 #ifndef ITK_MANUAL_INSTANTIATION
114 #include "rtkDrawGeometricPhantomImageFilter.hxx"
115 #endif
116 
117 #endif
Container for a geometric phantom, i.e., a set of ConvexShapes.
itk::Vector< ScalarType, Dimension > VectorType
itk::SmartPointer< const Self > ConstPointer
~DrawGeometricPhantomImageFilter() override=default
Draws a GeometricPhantom in a 3D image.
itk::Matrix< ScalarType, Dimension, Dimension > RotationMatrixType
#define itkSetMacro(name, type)