ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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 : m_X(p)
42 , m_Sign(s)
43 {}
44
45 Point1D(const Point1D & point) = default;
46
47 [[nodiscard]] double
48 getX() const
49 {
50 return m_X;
51 }
52
53 [[nodiscard]] int
54 getSign() const
55 {
56 return m_Sign;
57 }
58};
59
67template <typename TInputMesh, typename TOutputImage>
68class ITK_TEMPLATE_EXPORT TriangleMeshToBinaryImageFilter : public ImageSource<TOutputImage>
69{
70public:
71 ITK_DISALLOW_COPY_AND_MOVE(TriangleMeshToBinaryImageFilter);
72
78
79 using IndexType = typename TOutputImage::IndexType;
80 using SizeType = typename TOutputImage::SizeType;
81 using OutputImageType = TOutputImage;
82 using OutputImagePointer = typename OutputImageType::Pointer;
83 using ValueType = typename OutputImageType::ValueType;
84 using SpacingType = typename OutputImageType::SpacingType;
85 using DirectionType = typename OutputImageType::DirectionType;
86
88 itkNewMacro(Self);
89
91 itkOverrideGetNameOfClassMacro(TriangleMeshToBinaryImageFilter);
92
95
97 using InputMeshType = TInputMesh;
98 using InputMeshPointer = typename InputMeshType::Pointer;
99 using InputPointType = typename InputMeshType::PointType;
100 using InputPixelType = typename InputMeshType::PixelType;
101 using InputCellTraitsType = typename InputMeshType::MeshTraits::CellTraits;
102 using CellType = typename InputMeshType::CellType;
103 using CellsContainerPointer = typename InputMeshType::CellsContainerPointer;
104 using CellsContainerIterator = typename InputMeshType::CellsContainerIterator;
105
106 using InputPointsContainer = typename InputMeshType::PointsContainer;
107 using InputPointsContainerPointer = typename InputPointsContainer::Pointer;
108 using InputPointsContainerIterator = typename InputPointsContainer::Iterator;
109
112
115
117
118 using Point1DVector = std::vector<Point1D>;
119 using Point1DArray = std::vector<std::vector<Point1D>>;
120
121 using Point2DVector = std::vector<Point2DType>;
122 using Point2DArray = std::vector<std::vector<Point2DType>>;
123
124 using PointVector = std::vector<PointType>;
125 using PointArray = std::vector<std::vector<PointType>>;
126
132 itkSetMacro(Spacing, SpacingType);
133 virtual void
134 SetSpacing(const double spacing[3]);
136 virtual void
137 SetSpacing(const float spacing[3]);
138
139 itkGetConstReferenceMacro(Spacing, SpacingType);
140
145 itkSetMacro(Direction, DirectionType);
146 itkGetConstMacro(Direction, DirectionType);
154 itkSetMacro(InsideValue, ValueType);
155 itkGetConstMacro(InsideValue, ValueType);
163 itkSetMacro(OutsideValue, ValueType);
164 itkGetConstMacro(OutsideValue, ValueType);
171 itkSetMacro(Origin, PointType);
172 virtual void
173 SetOrigin(const double origin[3]);
175 virtual void
176 SetOrigin(const float origin[3]);
177
178 itkGetConstReferenceMacro(Origin, PointType);
179
182 itkSetMacro(Index, IndexType);
183 itkGetConstMacro(Index, IndexType);
187 itkSetMacro(Size, SizeType);
188 itkGetConstMacro(Size, SizeType);
192 void
194
195 void
197 {
198 if (InfoImage != m_InfoImage)
199 {
200 this->Modified();
201 m_InfoImage = InfoImage;
202 }
203 }
204
206 InputMeshType *
208
210 GetInput(unsigned int idx);
211
212 /* Set the tolerance for doing spatial searches of the polydata. */
213 itkSetMacro(Tolerance, double);
214 itkGetConstMacro(Tolerance, double);
215
216protected:
219
220 void
222 {} // do nothing
223 void
224 GenerateData() override;
225
226 virtual void
228
230 static int
231 PolygonToImageRaster(PointVector coords, Point1DArray & zymatrix, int extent[6]);
232
234
236
238
240
241 PointType m_Origin{}; // start value
242
243 double m_Tolerance{};
244
247
249
250 void
251 PrintSelf(std::ostream & os, Indent indent) const override;
252
253private:
254 static bool
256
257 static bool
259};
260} // end namespace itk
261
262#ifndef ITK_MANUAL_INSTANTIATION
263# include "itkTriangleMeshToBinaryImageFilter.hxx"
264#endif
265
266#endif
Array class with size defined at construction time.
Definition itkArray.h:48
typename OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
virtual void Modified() const
Point1D(const double p, const int s)
Point1D()=default
Point1D(const Point1D &point)=default
TPointsContainer PointsContainer
A superclass of the N-dimensional mesh structure; supports point (geometric coordinate and attribute)...
Definition itkPointSet.h:82
A templated class holding a geometric point in n-Dimensional space.
Definition itkPoint.h:54
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
Implements transparent reference counting.
typename InputPointsContainer::Iterator InputPointsContainerIterator
std::vector< std::vector< Point2DType > > Point2DArray
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
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
std::vector< std::vector< PointType > > PointArray
void SetInput(InputMeshType *input)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
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