ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkTriangleMeshToBinaryImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software 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 #ifndef itkTriangleMeshToBinaryImageFilter_h
19 #define itkTriangleMeshToBinaryImageFilter_h
20 
21 #include "itkImageSource.h"
22 
23 #include "itkPolygonCell.h"
24 #include "itkMapContainer.h"
25 #include "itkVectorContainer.h"
27 #include "itkPointSet.h"
28 
29 #include <vector>
30 
31 namespace itk
32 {
33 class Point1D
34 {
35 public:
36  double m_X;
37  int m_Sign;
38 
39  Point1D(){}
40  Point1D(const double p, const int s)
41  {
42  m_X = p;
43  m_Sign = s;
44  }
45 
46  Point1D(const Point1D & point)
47  {
48  m_X = point.m_X;
49  m_Sign = point.m_Sign;
50  }
51 
52  double getX() const
53  {
54  return m_X;
55  }
56 
57  int getSign() const
58  {
59  return m_Sign;
60  }
61 };
62 
70 template< typename TInputMesh, typename TOutputImage >
71 class TriangleMeshToBinaryImageFilter:public ImageSource< TOutputImage >
72 {
73 public:
79 
80  typedef typename TOutputImage::IndexType IndexType;
81  typedef typename TOutputImage::SizeType SizeType;
82  typedef TOutputImage OutputImageType;
83  typedef typename OutputImageType::Pointer OutputImagePointer;
84  typedef typename OutputImageType::ValueType ValueType;
85  typedef typename OutputImageType::SpacingType SpacingType;
86  typedef typename OutputImageType::DirectionType DirectionType;
87 
89  itkNewMacro(Self);
90 
93 
96 
98  typedef TInputMesh InputMeshType;
99  typedef typename InputMeshType::Pointer InputMeshPointer;
100  typedef typename InputMeshType::PointType InputPointType;
101  typedef typename InputMeshType::PixelType InputPixelType;
102  typedef typename InputMeshType::MeshTraits::CellTraits InputCellTraitsType;
103  typedef typename InputMeshType::CellType CellType;
104  typedef typename InputMeshType::CellsContainerPointer CellsContainerPointer;
105  typedef typename InputMeshType::CellsContainerIterator CellsContainerIterator;
106 
107  typedef typename InputMeshType::PointsContainer InputPointsContainer;
108  typedef typename InputPointsContainer::Pointer InputPointsContainerPointer;
109  typedef typename InputPointsContainer::Iterator InputPointsContainerIterator;
110 
113 
116 
118 
119  typedef std::vector< Point1D > Point1DVector;
120  typedef std::vector< std::vector< Point1D > > Point1DArray;
121 
122  typedef std::vector< Point2DType > Point2DVector;
123  typedef std::vector< std::vector< Point2DType > > Point2DArray;
124 
125  typedef std::vector< PointType > PointVector;
126  typedef std::vector< std::vector< PointType > > PointArray;
127 
128  typedef std::vector< int > StencilIndexVector;
133  itkSetMacro(Spacing, SpacingType);
134  virtual void SetSpacing(const double spacing[3]);
136 
137  virtual void SetSpacing(const float spacing[3]);
138 
139  itkGetConstReferenceMacro(Spacing, SpacingType);
140 
144  itkSetMacro(Direction, DirectionType);
145  itkGetConstMacro(Direction, DirectionType);
147 
153  itkSetMacro(InsideValue, ValueType);
154  itkGetConstMacro(InsideValue, ValueType);
156 
162  itkSetMacro(OutsideValue, ValueType);
163  itkGetConstMacro(OutsideValue, ValueType);
165 
170  itkSetMacro(Origin, PointType);
171  virtual void SetOrigin(const double origin[3]);
173 
174  virtual void SetOrigin(const float origin[3]);
175 
176  itkGetConstReferenceMacro(Origin, PointType);
177 
179  itkSetMacro(Index, IndexType);
180  itkGetConstMacro(Index, IndexType);
182 
184  itkSetMacro(Size, SizeType);
185  itkGetConstMacro(Size, SizeType);
187 
189  using Superclass::SetInput;
190  void SetInput(InputMeshType *input);
191 
192  void SetInfoImage(OutputImageType *InfoImage)
193  {
194  if ( InfoImage != m_InfoImage )
195  {
196  this->Modified();
197  m_InfoImage = InfoImage;
198  }
199  }
200 
203 
204  InputMeshType * GetInput(unsigned int idx);
205 
206  /* Set the tolerance for doing spatial searches of the polydata. */
207  itkSetMacro(Tolerance, double);
208  itkGetConstMacro(Tolerance, double);
209 
210 protected:
213 
214  virtual void GenerateOutputInformation() ITK_OVERRIDE {} // do nothing
215  virtual void GenerateData() ITK_OVERRIDE;
216 
217  virtual void RasterizeTriangles();
218 
219  static int PolygonToImageRaster(PointVector coords, Point1DArray & zymatrix, int extent[6]);
220 
222 
224 
226 
228 
229  PointType m_Origin; //start value
230 
231  double m_Tolerance;
232 
235 
237 
239 
240  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
241 
242 private:
243  TriangleMeshToBinaryImageFilter(const Self &) ITK_DELETE_FUNCTION;
244  void operator=(const Self &) ITK_DELETE_FUNCTION;
245 
246  static bool ComparePoints2D(Point2DType a, Point2DType b);
247 
248  static bool ComparePoints1D(Point1D a, Point1D b);
249 };
250 } // end namespace itk
251 
252 #ifndef ITK_MANUAL_INSTANTIATION
253 #include "itkTriangleMeshToBinaryImageFilter.hxx"
254 #endif
255 
256 #endif
Point1D(const Point1D &point)
InputMeshType::CellsContainerIterator CellsContainerIterator
3D Rasterization algorithm Courtesy of Dr David Gobbi of Atamai Inc.
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
void SetInput(InputMeshType *input)
Superclass::OutputImageRegionType OutputImageRegionType
virtual void SetOrigin(PointType _arg)
InputPointsContainer::Iterator InputPointsContainerIterator
InputMeshType::CellsContainerPointer CellsContainerPointer
Base class for all process objects that output image data.
Point1D(const double p, const int s)
virtual void SetSpacing(SpacingType _arg)
virtual void GenerateData() override
InputMeshType::MeshTraits::CellTraits InputCellTraitsType
static bool ComparePoints2D(Point2DType a, Point2DType b)
static int PolygonToImageRaster(PointVector coords, Point1DArray &zymatrix, int extent[6])
std::vector< std::vector< Point1D > > Point1DArray
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:84
MeshTraits::PointsContainer PointsContainer
Definition: itkPointSet.h:107
virtual void Modified() const
static bool ComparePoints1D(Point1D a, Point1D b)
std::vector< std::vector< PointType > > PointArray
OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
std::vector< std::vector< Point2DType > > Point2DArray
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:71
virtual void PrintSelf(std::ostream &os, Indent indent) const override