RTK  2.5.0
Reconstruction Toolkit
rtkForbildPhantomFileReader.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 rtkForbildPhantomFileReader_h
20 #define rtkForbildPhantomFileReader_h
21 
22 #include "RTKExport.h"
23 #include <itkLightProcessObject.h>
24 #include "rtkGeometricPhantom.h"
25 
26 namespace rtk
27 {
28 
43 {
44 public:
45  ITK_DISALLOW_COPY_AND_MOVE(ForbildPhantomFileReader);
46 
52 
54  static constexpr unsigned int Dimension = ConvexShape::Dimension;
61 
63  itkNewMacro(Self);
64 
66 #ifdef itkOverrideGetNameOfClassMacro
67  itkOverrideGetNameOfClassMacro(ForbildPhantomFileReader);
68 #else
70 #endif
71 
72 
74  itkGetModifiableObjectMacro(GeometricPhantom, GeometricPhantom);
75  itkSetObjectMacro(GeometricPhantom, GeometricPhantom);
77 
79  itkGetStringMacro(Filename);
80  itkSetStringMacro(Filename);
82 
84  virtual void
85  GenerateOutputInformation();
86 
87 protected:
89  ForbildPhantomFileReader() = default;
90 
92  ~ForbildPhantomFileReader() override = default;
93 
94  void
95  CreateForbildSphere(const std::string & s);
96  void
97  CreateForbildBox(const std::string & s);
98  void
99  CreateForbildCylinder(const std::string & s, const std::string & fig);
100  void
101  CreateForbildElliptCyl(const std::string & s, const std::string & fig);
102  void
103  CreateForbildEllipsoid(const std::string & s, const std::string & fig);
104  void
105  CreateForbildCone(const std::string & s, const std::string & fig);
106  void
107  CreateForbildTetrahedron(const std::string & s);
109  ComputeRotationMatrixBetweenVectors(const VectorType & source, const VectorType & dest) const;
110 
111  bool
112  FindParameterInString(const std::string & name, const std::string & s, ScalarType & param);
113  bool
114  FindVectorInString(const std::string & name, const std::string & s, VectorType & vec);
115  void
116  FindClipPlanes(const std::string & s);
117  void
118  FindUnions(const std::string & s);
119 
120 private:
122  std::string m_Filename;
126  std::vector<int> m_UnionWith;
127 };
128 
129 } // end namespace rtk
130 
131 #endif
Container for a geometric phantom, i.e., a set of ConvexShapes.
itk::Vector< ScalarType, Dimension > VectorType
static constexpr unsigned int Dimension
constexpr unsigned int Dimension
itk::SmartPointer< Self > Pointer
std::vector< ConvexShapePointer > ConvexShapeVector
GeometricPhantom::ConvexShapeVector ConvexShapeVectorType
Reads GeometricPhantom from a Forbild file.
itk::Vector< ScalarType, Dimension > PointType
itk::Matrix< ScalarType, Dimension, Dimension > RotationMatrixType