ITK  6.0.0
Insight Toolkit
itkBinaryMask3DMeshSource.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 itkBinaryMask3DMeshSource_h
19#define itkBinaryMask3DMeshSource_h
20
21#include "vnl/vnl_matrix_fixed.h"
22#include "itkMesh.h"
24#include "itkTriangleCell.h"
25#include "itkCovariantVector.h"
28
29namespace itk
30{
70template <typename TInputImage, typename TOutputMesh>
71class ITK_TEMPLATE_EXPORT BinaryMask3DMeshSource : public ImageToMeshFilter<TInputImage, TOutputMesh>
72{
73public:
74 ITK_DISALLOW_COPY_AND_MOVE(BinaryMask3DMeshSource);
75
81
83 itkNewMacro(Self);
84
86 itkOverrideGetNameOfClassMacro(BinaryMask3DMeshSource);
87
89 using OutputMeshType = TOutputMesh;
90 using OMeshTraits = typename OutputMeshType::MeshTraits;
92 using OPixelType = typename OMeshTraits::PixelType;
93
96 using CellTraits = typename OutputMeshType::CellTraits;
97 using PointsContainerPointer = typename OutputMeshType::PointsContainerPointer;
98 using PointsContainer = typename OutputMeshType::PointsContainer;
99 using CellsContainerPointer = typename OutputMeshType::CellsContainerPointer;
100 using CellsContainer = typename OutputMeshType::CellsContainer;
103
108 using TriCellAutoPointer = typename TriCell::SelfAutoPointer;
109
111 using InputImageType = TInputImage;
114 using InputPixelType = typename InputImageType::PixelType;
115 using SpacingType = typename InputImageType::SpacingType;
119
122
124
127
128 itkSetMacro(ObjectValue, InputPixelType);
129
130 itkGetConstMacro(NumberOfNodes, SizeValueType);
131 itkGetConstMacro(NumberOfCells, SizeValueType);
132
134 using Superclass::SetInput;
135 virtual void
136 SetInput(const InputImageType * image);
137
138 void
140 {
141 if (iRegion != m_RegionOfInterest)
142 {
143 this->m_RegionOfInterest = iRegion;
144 this->m_RegionOfInterestProvidedByUser = true;
145 this->Modified();
146 }
147 }
148
149 itkGetConstReferenceMacro(RegionOfInterest, RegionType);
150
151protected:
154 void
155 PrintSelf(std::ostream & os, Indent indent) const override;
156
157 void
158 GenerateData() override;
159
160
161 bool m_RegionOfInterestProvidedByUser{ false };
162 RegionType m_RegionOfInterest{};
163
164 void
166 {} // do nothing override
167
168private:
170
171 void
173
174 void
175 XFlip(unsigned char * x); // 7 kinds of transformation
176
177 void
178 YFlip(unsigned char * x);
179
180 void
181 ZFlip(unsigned char * x);
182
183 void
184 XRotation(unsigned char * x);
185
186 void
187 YRotation(unsigned char * x);
188
189 void
190 ZRotation(unsigned char * x);
191
192 void
193 inverse(unsigned char * x);
194
195 void
196 InitializeLUT(); // initialize the look up table before the mesh
197 // construction
198
199 void
200 AddCells(unsigned char celltype, unsigned char celltran, int index);
201
202 void
203 AddNodes(int index,
204 unsigned char * nodesid,
205 IdentifierType * globalnodesid,
206 IdentifierType ** currentrowtmp,
207 IdentifierType ** currentframetmp);
208
209 void
210 CellTransfer(unsigned char * nodesid, unsigned char celltran);
211
213 SearchThroughLastRow(int index, int start, int end);
214
216 SearchThroughLastFrame(int index, int start, int end);
217
218 unsigned char m_LUT[256][2]{}; // the two lookup tables
219
220 IdentifierType m_LastVoxel[14]{};
221 IdentifierType m_CurrentVoxel[14]{};
222
223 IdentifierType ** m_LastRow{ nullptr };
224 IdentifierType ** m_LastFrame{ nullptr };
225 IdentifierType ** m_CurrentRow{ nullptr };
226 IdentifierType ** m_CurrentFrame{ nullptr };
227
228 int m_CurrentRowIndex{ 0 };
229 int m_CurrentFrameIndex{ 0 };
230 int m_LastRowNum{ 0 };
231 int m_LastFrameNum{ 0 };
232 int m_CurrentRowNum{ 200 };
233 int m_CurrentFrameNum{ 2000 };
234 unsigned char m_AvailableNodes[14]{};
235
236 double m_LocationOffset[14][3]{};
237
238 SizeValueType m_NumberOfNodes{ 0 };
239 SizeValueType m_NumberOfCells{ 0 };
240
241 int m_NodeLimit{ 2000 };
242 int m_CellLimit{ 4000 };
243 int m_ImageWidth{ 0 };
244 int m_ImageHeight{ 0 };
245 int m_ImageDepth{ 0 };
246 int m_ColFlag{ 0 };
247 int m_RowFlag{ 0 };
248 int m_FrameFlag{ 0 };
249 int m_LastRowIndex{ 0 };
250 int m_LastVoxelIndex{ 0 };
251 int m_LastFrameIndex{ 0 };
252
253 unsigned char m_PointFound{ 0 };
254 InputPixelType m_ObjectValue{};
255
259 OutputMeshType * m_OutputMesh{};
260 const InputImageType * m_InputImage{};
261};
262} // end namespace itk
263
264#ifndef ITK_MANUAL_INSTANTIATION
265# include "itkBinaryMask3DMeshSource.hxx"
266#endif
267
268#endif
typename OutputMeshType::PointsContainerPointer PointsContainerPointer
typename InputImageType::Pointer InputImagePointer
typename OutputMeshType::PointType OPointType
typename InputImageType::RegionType RegionType
typename InputImageType::PointType OriginType
typename OutputMeshType::CellsContainerPointer CellsContainerPointer
typename OutputMeshType::PointsContainer PointsContainer
void XRotation(unsigned char *x)
virtual void SetInput(const InputImageType *image)
void YRotation(unsigned char *x)
void inverse(unsigned char *x)
typename OutputMeshType::MeshTraits OMeshTraits
typename InputImageType::ConstPointer InputImageConstPointer
void SetRegionOfInterest(const RegionType &iRegion)
typename OutputMeshType::CellTraits CellTraits
void XFlip(unsigned char *x)
void ZRotation(unsigned char *x)
typename InputImageType::IndexType InputImageIndexType
void GenerateData() override
typename InputImageType::SizeType SizeType
typename OutputMeshType::CellsContainer CellsContainer
typename TriCell::SelfAutoPointer TriCellAutoPointer
typename InputImageType::SpacingType SpacingType
typename InputImageType::SizeType InputImageSizeType
IdentifierType SearchThroughLastFrame(int index, int start, int end)
typename OMeshTraits::PixelType OPixelType
typename OutputMeshType::Pointer OutputMeshPointer
typename InputImageType::PixelType InputPixelType
IdentifierType SearchThroughLastRow(int index, int start, int end)
void CellTransfer(unsigned char *nodesid, unsigned char celltran)
void PrintSelf(std::ostream &os, Indent indent) const override
void YFlip(unsigned char *x)
void AddNodes(int index, unsigned char *nodesid, IdentifierType *globalnodesid, IdentifierType **currentrowtmp, IdentifierType **currentframetmp)
void AddCells(unsigned char celltype, unsigned char celltran, int index)
void ZFlip(unsigned char *x)
An abstract interface for cells.
A templated class holding a n-Dimensional covariant vector.
A multi-dimensional iterator templated over image type that walks a region of pixels.
ImageToMeshFilter is the base class for all process objects that output Mesh data and require image d...
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
SizeValueType IdentifierType
Definition: itkIntTypes.h:90
unsigned long SizeValueType
Definition: itkIntTypes.h:86