ITK  6.0.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 * 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#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
31namespace itk
32{
34{
35public:
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
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
72template <typename TInputMesh, typename TOutputImage>
73class ITK_TEMPLATE_EXPORT TriangleMeshToBinaryImageFilter : public ImageSource<TOutputImage>
74{
75public:
76 ITK_DISALLOW_COPY_AND_MOVE(TriangleMeshToBinaryImageFilter);
77
83
86 using OutputImageType = TOutputImage;
88 using ValueType = typename OutputImageType::ValueType;
89 using SpacingType = typename OutputImageType::SpacingType;
91
93 itkNewMacro(Self);
94
96 itkOverrideGetNameOfClassMacro(TriangleMeshToBinaryImageFilter);
97
99 using typename Superclass::OutputImageRegionType;
100
102 using InputMeshType = TInputMesh;
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;
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
136 itkSetMacro(Spacing, SpacingType);
137 virtual void
138 SetSpacing(const double spacing[3]);
141 virtual void
142 SetSpacing(const float spacing[3]);
143
144 itkGetConstReferenceMacro(Spacing, SpacingType);
145
149 itkSetMacro(Direction, DirectionType);
150 itkGetConstMacro(Direction, DirectionType);
158 itkSetMacro(InsideValue, ValueType);
159 itkGetConstMacro(InsideValue, ValueType);
167 itkSetMacro(OutsideValue, ValueType);
168 itkGetConstMacro(OutsideValue, ValueType);
175 itkSetMacro(Origin, PointType);
176 virtual void
177 SetOrigin(const double origin[3]);
180 virtual void
181 SetOrigin(const float origin[3]);
182
183 itkGetConstReferenceMacro(Origin, PointType);
184
186 itkSetMacro(Index, IndexType);
187 itkGetConstMacro(Index, IndexType);
191 itkSetMacro(Size, SizeType);
192 itkGetConstMacro(Size, SizeType);
196 using Superclass::SetInput;
197 void
199
200 void
202 {
203 if (InfoImage != m_InfoImage)
204 {
205 this->Modified();
206 m_InfoImage = InfoImage;
207 }
208 }
209
211 InputMeshType *
213
215 GetInput(unsigned int idx);
216
217 /* Set the tolerance for doing spatial searches of the polydata. */
218 itkSetMacro(Tolerance, double);
219 itkGetConstMacro(Tolerance, double);
220
221protected:
224
225 void
227 {} // do nothing
228 void
229 GenerateData() override;
230
231 virtual void
233
235 static int
236 PolygonToImageRaster(PointVector coords, Point1DArray & zymatrix, int extent[6]);
237
238 OutputImageType * m_InfoImage{};
239
240 IndexType m_Index{};
241
242 SizeType m_Size{};
243
244 SpacingType m_Spacing{};
245
246 PointType m_Origin{}; // start value
247
248 double m_Tolerance{};
249
250 ValueType m_InsideValue{};
251 ValueType m_OutsideValue{};
252
253 DirectionType m_Direction{};
254
255 void
256 PrintSelf(std::ostream & os, Indent indent) const override;
257
258private:
259 static bool
261
262 static bool
264};
265} // end namespace itk
266
267#ifndef ITK_MANUAL_INSTANTIATION
268# include "itkTriangleMeshToBinaryImageFilter.hxx"
269#endif
270
271#endif
Base class for all process objects that output image data.
typename OutputImageType::Pointer OutputImagePointer
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Point1D(const double p, const int s)
Point1D()=default
Point1D(const Point1D &point)
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition: itkPointSet.h:82
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
3D Rasterization algorithm Courtesy of Dr David Gobbi of Atamai Inc.
typename InputPointsContainer::Iterator InputPointsContainerIterator
typename InputMeshType::PointsContainer InputPointsContainer
typename OutputImageType::DirectionType DirectionType
typename OutputImageType::SpacingType SpacingType
virtual void SetOrigin(const double origin[3])
InputMeshType * GetInput(unsigned int idx)
virtual void SetSpacing(const float spacing[3])
void PrintSelf(std::ostream &os, Indent indent) const override
typename InputMeshType::CellsContainerPointer CellsContainerPointer
typename InputPointsContainer::Pointer InputPointsContainerPointer
static bool ComparePoints2D(Point2DType a, Point2DType b)
virtual void SetOrigin(const float origin[3])
~TriangleMeshToBinaryImageFilter() override=default
std::vector< std::vector< PointType > > PointArray
typename InputMeshType::MeshTraits::CellTraits InputCellTraitsType
typename PointSetType::PointsContainer PointsContainer
static int PolygonToImageRaster(PointVector coords, Point1DArray &zymatrix, int extent[6])
typename InputMeshType::CellsContainerIterator CellsContainerIterator
virtual void SetSpacing(const double spacing[3])
static bool ComparePoints1D(Point1D a, Point1D b)
std::vector< std::vector< Point1D > > Point1DArray
void SetInput(InputMeshType *input)
std::vector< std::vector< Point2DType > > Point2DArray
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
*par Constraints *The filter image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
Represent a n-dimensional index in a n-dimensional image.
Definition: itkIndex.h:69
Represent a n-dimensional size (bounds) of a n-dimensional image.
Definition: itkSize.h:70