ITK  5.2.0
Insight Toolkit
itkTriangleMeshToBinaryImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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() = default;
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
53  getX() const
54  {
55  return m_X;
56  }
57 
58  int
59  getSign() const
60  {
61  return m_Sign;
62  }
63 };
64 
72 template <typename TInputMesh, typename TOutputImage>
73 class ITK_TEMPLATE_EXPORT TriangleMeshToBinaryImageFilter : public ImageSource<TOutputImage>
74 {
75 public:
76  ITK_DISALLOW_COPY_AND_MOVE(TriangleMeshToBinaryImageFilter);
77 
83 
85  using SizeType = typename TOutputImage::SizeType;
86  using OutputImageType = TOutputImage;
87  using OutputImagePointer = typename OutputImageType::Pointer;
88  using ValueType = typename OutputImageType::ValueType;
89  using SpacingType = typename OutputImageType::SpacingType;
91 
93  itkNewMacro(Self);
94 
97 
99  using OutputImageRegionType = typename Superclass::OutputImageRegionType;
100 
102  using InputMeshType = TInputMesh;
103  using InputMeshPointer = typename InputMeshType::Pointer;
105  using InputPixelType = typename InputMeshType::PixelType;
106  using InputCellTraitsType = typename InputMeshType::MeshTraits::CellTraits;
107  using CellType = typename InputMeshType::CellType;
108  using CellsContainerPointer = typename InputMeshType::CellsContainerPointer;
109  using CellsContainerIterator = typename InputMeshType::CellsContainerIterator;
110 
111  using InputPointsContainer = typename InputMeshType::PointsContainer;
112  using InputPointsContainerPointer = typename InputPointsContainer::Pointer;
113  using InputPointsContainerIterator = typename InputPointsContainer::Iterator;
114 
117 
120 
122 
123  using Point1DVector = std::vector<Point1D>;
124  using Point1DArray = std::vector<std::vector<Point1D>>;
125 
126  using Point2DVector = std::vector<Point2DType>;
127  using Point2DArray = std::vector<std::vector<Point2DType>>;
128 
129  using PointVector = std::vector<PointType>;
130  using PointArray = std::vector<std::vector<PointType>>;
131 
132  using StencilIndexVector = std::vector<int>;
137  itkSetMacro(Spacing, SpacingType);
138  virtual void
139  SetSpacing(const double spacing[3]);
141 
142  virtual void
143  SetSpacing(const float spacing[3]);
144 
145  itkGetConstReferenceMacro(Spacing, SpacingType);
146 
150  itkSetMacro(Direction, DirectionType);
151  itkGetConstMacro(Direction, DirectionType);
153 
159  itkSetMacro(InsideValue, ValueType);
160  itkGetConstMacro(InsideValue, ValueType);
162 
168  itkSetMacro(OutsideValue, ValueType);
169  itkGetConstMacro(OutsideValue, ValueType);
171 
176  itkSetMacro(Origin, PointType);
177  virtual void
178  SetOrigin(const double origin[3]);
180 
181  virtual void
182  SetOrigin(const float origin[3]);
183 
184  itkGetConstReferenceMacro(Origin, PointType);
185 
187  itkSetMacro(Index, IndexType);
188  itkGetConstMacro(Index, IndexType);
190 
192  itkSetMacro(Size, SizeType);
193  itkGetConstMacro(Size, SizeType);
195 
197  using Superclass::SetInput;
198  void
199  SetInput(InputMeshType * input);
200 
201  void
203  {
204  if (InfoImage != m_InfoImage)
205  {
206  this->Modified();
207  m_InfoImage = InfoImage;
208  }
209  }
210 
212  InputMeshType *
213  GetInput();
214 
215  InputMeshType *
216  GetInput(unsigned int idx);
217 
218  /* Set the tolerance for doing spatial searches of the polydata. */
219  itkSetMacro(Tolerance, double);
220  itkGetConstMacro(Tolerance, double);
221 
222 protected:
224  ~TriangleMeshToBinaryImageFilter() override = default;
225 
226  void
228  {} // do nothing
229  void
230  GenerateData() override;
231 
232  virtual void
233  RasterizeTriangles();
234 
235  static int
236  PolygonToImageRaster(PointVector coords, Point1DArray & zymatrix, int extent[6]);
237 
239 
241 
243 
245 
246  PointType m_Origin; // start value
247 
248  double m_Tolerance;
249 
252 
254 
256 
257  void
258  PrintSelf(std::ostream & os, Indent indent) const override;
259 
260 private:
261  static bool
262  ComparePoints2D(Point2DType a, Point2DType b);
263 
264  static bool
265  ComparePoints1D(Point1D a, Point1D b);
266 };
267 } // end namespace itk
268 
269 #ifndef ITK_MANUAL_INSTANTIATION
270 # include "itkTriangleMeshToBinaryImageFilter.hxx"
271 #endif
272 
273 #endif
itk::TriangleMeshToBinaryImageFilter::m_InfoImage
OutputImageType * m_InfoImage
Definition: itkTriangleMeshToBinaryImageFilter.h:238
itk::TriangleMeshToBinaryImageFilter::DirectionType
typename OutputImageType::DirectionType DirectionType
Definition: itkTriangleMeshToBinaryImageFilter.h:90
itk::Index
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:66
itk::TriangleMeshToBinaryImageFilter::CellType
typename InputMeshType::CellType CellType
Definition: itkTriangleMeshToBinaryImageFilter.h:107
itk::TriangleMeshToBinaryImageFilter::CellsContainerIterator
typename InputMeshType::CellsContainerIterator CellsContainerIterator
Definition: itkTriangleMeshToBinaryImageFilter.h:109
itk::GTest::TypedefsAndConstructors::Dimension2::DirectionType
ImageBaseType::DirectionType DirectionType
Definition: itkGTestTypedefsAndConstructors.h:52
itk::TriangleMeshToBinaryImageFilter::CellsContainerPointer
typename InputMeshType::CellsContainerPointer CellsContainerPointer
Definition: itkTriangleMeshToBinaryImageFilter.h:108
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itk::PointSet
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:82
itkImageSource.h
itk::TriangleMeshToBinaryImageFilter::m_Tolerance
double m_Tolerance
Definition: itkTriangleMeshToBinaryImageFilter.h:248
itk::TriangleMeshToBinaryImageFilter::m_Index
IndexType m_Index
Definition: itkTriangleMeshToBinaryImageFilter.h:240
itk::Size
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:69
itk::TriangleMeshToBinaryImageFilter::InputPixelType
typename InputMeshType::PixelType InputPixelType
Definition: itkTriangleMeshToBinaryImageFilter.h:105
itk::TriangleMeshToBinaryImageFilter::m_OutsideValue
ValueType m_OutsideValue
Definition: itkTriangleMeshToBinaryImageFilter.h:251
itk::Point1D
Definition: itkTriangleMeshToBinaryImageFilter.h:33
itk::TriangleMeshToBinaryImageFilter::StencilIndexVector
std::vector< int > StencilIndexVector
Definition: itkTriangleMeshToBinaryImageFilter.h:132
itk::TriangleMeshToBinaryImageFilter::InputPointsContainer
typename InputMeshType::PointsContainer InputPointsContainer
Definition: itkTriangleMeshToBinaryImageFilter.h:111
itk::TriangleMeshToBinaryImageFilter::InputPointsContainerPointer
typename InputPointsContainer::Pointer InputPointsContainerPointer
Definition: itkTriangleMeshToBinaryImageFilter.h:112
itk::TriangleMeshToBinaryImageFilter::m_Spacing
SpacingType m_Spacing
Definition: itkTriangleMeshToBinaryImageFilter.h:244
itk::TriangleMeshToBinaryImageFilter::InputMeshPointer
typename InputMeshType::Pointer InputMeshPointer
Definition: itkTriangleMeshToBinaryImageFilter.h:103
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::TriangleMeshToBinaryImageFilter::InputPointsContainerIterator
typename InputPointsContainer::Iterator InputPointsContainerIterator
Definition: itkTriangleMeshToBinaryImageFilter.h:113
itk::TriangleMeshToBinaryImageFilter::SizeType
typename TOutputImage::SizeType SizeType
Definition: itkTriangleMeshToBinaryImageFilter.h:85
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::SmartPointer< Self >
itk::Point1D::getSign
int getSign() const
Definition: itkTriangleMeshToBinaryImageFilter.h:59
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::TriangleMeshToBinaryImageFilter::SetInfoImage
void SetInfoImage(OutputImageType *InfoImage)
Definition: itkTriangleMeshToBinaryImageFilter.h:202
itk::TriangleMeshToBinaryImageFilter::m_Size
SizeType m_Size
Definition: itkTriangleMeshToBinaryImageFilter.h:242
itk::TriangleMeshToBinaryImageFilter::SpacingType
typename OutputImageType::SpacingType SpacingType
Definition: itkTriangleMeshToBinaryImageFilter.h:89
itk::TriangleMeshToBinaryImageFilter::PointArray
std::vector< std::vector< PointType > > PointArray
Definition: itkTriangleMeshToBinaryImageFilter.h:130
itk::Point1D::Point1D
Point1D(const double p, const int s)
Definition: itkTriangleMeshToBinaryImageFilter.h:40
itk::TriangleMeshToBinaryImageFilter::m_Origin
PointType m_Origin
Definition: itkTriangleMeshToBinaryImageFilter.h:246
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::Point1D::m_X
double m_X
Definition: itkTriangleMeshToBinaryImageFilter.h:36
itk::TriangleMeshToBinaryImageFilter::Point2DVector
std::vector< Point2DType > Point2DVector
Definition: itkTriangleMeshToBinaryImageFilter.h:126
itk::TriangleMeshToBinaryImageFilter
3D Rasterization algorithm Courtesy of Dr David Gobbi of Atamai Inc.
Definition: itkTriangleMeshToBinaryImageFilter.h:73
itk::TriangleMeshToBinaryImageFilter::PointsContainer
typename PointSetType::PointsContainer PointsContainer
Definition: itkTriangleMeshToBinaryImageFilter.h:116
itk::TriangleMeshToBinaryImageFilter::ValueType
typename OutputImageType::ValueType ValueType
Definition: itkTriangleMeshToBinaryImageFilter.h:88
itk::TriangleMeshToBinaryImageFilter::PointVector
std::vector< PointType > PointVector
Definition: itkTriangleMeshToBinaryImageFilter.h:129
itk::TriangleMeshToBinaryImageFilter::OutputImageType
TOutputImage OutputImageType
Definition: itkTriangleMeshToBinaryImageFilter.h:86
itk::Point1D::m_Sign
int m_Sign
Definition: itkTriangleMeshToBinaryImageFilter.h:37
itk::TriangleMeshToBinaryImageFilter::m_Direction
DirectionType m_Direction
Definition: itkTriangleMeshToBinaryImageFilter.h:253
itkAutomaticTopologyMeshSource.h
itk::PointSet::PointsContainer
typename MeshTraits::PointsContainer PointsContainer
Definition: itkPointSet.h:107
itk::ImageSource::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkImageSource.h:92
itkVectorContainer.h
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itk::TriangleMeshToBinaryImageFilter::Point1DArray
std::vector< std::vector< Point1D > > Point1DArray
Definition: itkTriangleMeshToBinaryImageFilter.h:124
itkPolygonCell.h
itk::TriangleMeshToBinaryImageFilter::Point1DVector
std::vector< Point1D > Point1DVector
Definition: itkTriangleMeshToBinaryImageFilter.h:123
itk::Array< double >
itk::Point1D::Point1D
Point1D(const Point1D &point)
Definition: itkTriangleMeshToBinaryImageFilter.h:46
itk::TriangleMeshToBinaryImageFilter::m_InsideValue
ValueType m_InsideValue
Definition: itkTriangleMeshToBinaryImageFilter.h:250
itk::Point< double, 3 >
itk::TriangleMeshToBinaryImageFilter::Point2DArray
std::vector< std::vector< Point2DType > > Point2DArray
Definition: itkTriangleMeshToBinaryImageFilter.h:127
itk::Point1D::getX
double getX() const
Definition: itkTriangleMeshToBinaryImageFilter.h:53
itk::TriangleMeshToBinaryImageFilter::GenerateOutputInformation
void GenerateOutputInformation() override
Definition: itkTriangleMeshToBinaryImageFilter.h:227
itkPointSet.h
itk::TriangleMeshToBinaryImageFilter::IndexType
typename TOutputImage::IndexType IndexType
Definition: itkTriangleMeshToBinaryImageFilter.h:84
itk::TriangleMeshToBinaryImageFilter::m_StencilIndex
StencilIndexVector m_StencilIndex
Definition: itkTriangleMeshToBinaryImageFilter.h:255
itk::TriangleMeshToBinaryImageFilter::InputPointType
typename InputMeshType::PointType InputPointType
Definition: itkTriangleMeshToBinaryImageFilter.h:104
itkMapContainer.h
itk::TriangleMeshToBinaryImageFilter::InputCellTraitsType
typename InputMeshType::MeshTraits::CellTraits InputCellTraitsType
Definition: itkTriangleMeshToBinaryImageFilter.h:106
itk::Point1D::Point1D
Point1D()=default
itk::TriangleMeshToBinaryImageFilter::InputMeshType
TInputMesh InputMeshType
Definition: itkTriangleMeshToBinaryImageFilter.h:102