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"
23 #include "itkImageToMeshFilter.h"
24 #include "itkTriangleCell.h"
25 #include "itkCovariantVector.h"
28 
29 namespace itk
30 {
70 template <typename TInputImage, typename TOutputMesh>
71 class ITK_TEMPLATE_EXPORT BinaryMask3DMeshSource : public ImageToMeshFilter<TInputImage, TOutputMesh>
72 {
73 public:
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 
151 protected:
153  ~BinaryMask3DMeshSource() override;
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 
168 private:
170 
171  void
172  CreateMesh();
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
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::BinaryMask3DMeshSource::OMeshTraits
typename OutputMeshType::MeshTraits OMeshTraits
Definition: itkBinaryMask3DMeshSource.h:90
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::BinaryMask3DMeshSource::OPixelType
typename OMeshTraits::PixelType OPixelType
Definition: itkBinaryMask3DMeshSource.h:92
itk::BinaryMask3DMeshSource::InputImageSizeType
typename InputImageType::SizeType InputImageSizeType
Definition: itkBinaryMask3DMeshSource.h:169
itkCovariantVector.h
itk::BinaryMask3DMeshSource
Definition: itkBinaryMask3DMeshSource.h:71
itk::GTest::TypedefsAndConstructors::Dimension2::PointType
ImageBaseType::PointType PointType
Definition: itkGTestTypedefsAndConstructors.h:51
itk::BinaryMask3DMeshSource::GenerateOutputInformation
void GenerateOutputInformation() override
Definition: itkBinaryMask3DMeshSource.h:165
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::BinaryMask3DMeshSource::PointsContainerPointer
typename OutputMeshType::PointsContainerPointer PointsContainerPointer
Definition: itkBinaryMask3DMeshSource.h:97
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itkImageRegionConstIterator.h
itk::BinaryMask3DMeshSource::SizeValueType
itk::SizeValueType SizeValueType
Definition: itkBinaryMask3DMeshSource.h:126
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::BinaryMask3DMeshSource::InputImageIndexType
typename InputImageType::IndexType InputImageIndexType
Definition: itkBinaryMask3DMeshSource.h:121
itk::BinaryMask3DMeshSource::TriCellAutoPointer
typename TriCell::SelfAutoPointer TriCellAutoPointer
Definition: itkBinaryMask3DMeshSource.h:108
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::BinaryMask3DMeshSource::RegionType
typename InputImageType::RegionType RegionType
Definition: itkBinaryMask3DMeshSource.h:117
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itkMesh.h
itk::BinaryMask3DMeshSource::OriginType
typename InputImageType::PointType OriginType
Definition: itkBinaryMask3DMeshSource.h:116
itk::TriangleCell
Definition: itkTriangleCell.h:46
itk::BinaryMask3DMeshSource::SizeType
typename InputImageType::SizeType SizeType
Definition: itkBinaryMask3DMeshSource.h:118
itk::ImageToMeshFilter
ImageToMeshFilter is the base class for all process objects that output Mesh data and require image d...
Definition: itkImageToMeshFilter.h:36
itk::BinaryMask3DMeshSource::CellsContainer
typename OutputMeshType::CellsContainer CellsContainer
Definition: itkBinaryMask3DMeshSource.h:100
itk::BinaryMask3DMeshSource::InputImageConstPointer
typename InputImageType::ConstPointer InputImageConstPointer
Definition: itkBinaryMask3DMeshSource.h:113
itk::BinaryMask3DMeshSource::IdentifierType
itk::IdentifierType IdentifierType
Definition: itkBinaryMask3DMeshSource.h:125
itkImageToMeshFilter.h
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::BinaryMask3DMeshSource::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkBinaryMask3DMeshSource.h:112
itk::BinaryMask3DMeshSource::CellTraits
typename OutputMeshType::CellTraits CellTraits
Definition: itkBinaryMask3DMeshSource.h:96
itk::BinaryMask3DMeshSource::InputPixelType
typename InputImageType::PixelType InputPixelType
Definition: itkBinaryMask3DMeshSource.h:114
itk::CellInterface
An abstract interface for cells.
Definition: itkCellInterface.h:97
itk::BinaryMask3DMeshSource::SetRegionOfInterest
void SetRegionOfInterest(const RegionType &iRegion)
Definition: itkBinaryMask3DMeshSource.h:139
itk::BinaryMask3DMeshSource::OutputMeshPointer
typename OutputMeshType::Pointer OutputMeshPointer
Definition: itkBinaryMask3DMeshSource.h:95
itk::BinaryMask3DMeshSource::SpacingType
typename InputImageType::SpacingType SpacingType
Definition: itkBinaryMask3DMeshSource.h:115
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::BinaryMask3DMeshSource::CellsContainerPointer
typename OutputMeshType::CellsContainerPointer CellsContainerPointer
Definition: itkBinaryMask3DMeshSource.h:99
itk::BinaryMask3DMeshSource::InputImageType
TInputImage InputImageType
Definition: itkBinaryMask3DMeshSource.h:111
itkDefaultStaticMeshTraits.h
itk::ImageRegionConstIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionConstIterator.h:109
itk::BinaryMask3DMeshSource::OutputMeshType
TOutputMesh OutputMeshType
Definition: itkBinaryMask3DMeshSource.h:89
itkTriangleCell.h
itk::BinaryMask3DMeshSource::PointsContainer
typename OutputMeshType::PointsContainer PointsContainer
Definition: itkBinaryMask3DMeshSource.h:98
itk::IdentifierType
SizeValueType IdentifierType
Definition: itkIntTypes.h:90
itk::BinaryMask3DMeshSource::OPointType
typename OutputMeshType::PointType OPointType
Definition: itkBinaryMask3DMeshSource.h:91
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86