ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkBinaryMask3DMeshSource.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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 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_ASSIGN(BinaryMask3DMeshSource);
75 
81 
83  itkNewMacro(Self);
84 
87 
89  using OutputMeshType = TOutputMesh;
90  using OMeshTraits = typename OutputMeshType::MeshTraits;
92  using OPixelType = typename OMeshTraits::PixelType;
93 
95  using OutputMeshPointer = typename OutputMeshType::Pointer;
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;
112  using InputImagePointer = typename InputImageType::Pointer;
113  using InputImageConstPointer = typename InputImageType::ConstPointer;
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 SetInput(const InputImageType *inputImage);
136 
137  void SetRegionOfInterest( const RegionType & iRegion )
138  {
139  if( iRegion != m_RegionOfInterest )
140  {
141  this->m_RegionOfInterest = iRegion;
142  this->m_RegionOfInterestProvidedByUser = true;
143  this->Modified();
144  }
145  }
146 
147  itkGetConstReferenceMacro(RegionOfInterest, RegionType);
148 
149 protected:
151  ~BinaryMask3DMeshSource() override;
152  void PrintSelf(std::ostream & os, Indent indent) const override;
153 
154  void GenerateData() override;
155 
156 
157  bool m_RegionOfInterestProvidedByUser{false};
159 
160  void GenerateOutputInformation() override {} // do nothing override
161 
162 private:
164 
165  void CreateMesh();
166 
167  void XFlip(unsigned char *tp); // 7 kinds of transformation
168 
169  void YFlip(unsigned char *tp);
170 
171  void ZFlip(unsigned char *tp);
172 
173  void XRotation(unsigned char *tp);
174 
175  void YRotation(unsigned char *tp);
176 
177  void ZRotation(unsigned char *tp);
178 
179  void inverse(unsigned char *tp);
180 
181  void InitializeLUT(); // initialize the look up table before the mesh
182  // construction
183 
184  void AddCells(unsigned char celltype, unsigned char celltran, int index);
185 
186  void AddNodes(int index,
187  unsigned char *nodesid,
188  IdentifierType *globalnodesid,
189  IdentifierType **currentrowtmp,
190  IdentifierType **currentframetmp);
191 
192  void CellTransfer(unsigned char *nodesid, unsigned char celltran);
193 
194  IdentifierType SearchThroughLastRow(int index, int start, int end);
195 
196  IdentifierType SearchThroughLastFrame(int index, int start, int end);
197 
198  unsigned char m_LUT[256][2]; // the two lookup tables
199 
200  IdentifierType m_LastVoxel[14];
201  IdentifierType m_CurrentVoxel[14];
202 
203  IdentifierType **m_LastRow{nullptr};
204  IdentifierType **m_LastFrame{nullptr};
205  IdentifierType **m_CurrentRow{nullptr};
206  IdentifierType **m_CurrentFrame{nullptr};
207 
208  unsigned short m_CurrentRowIndex{0};
209  unsigned short m_CurrentFrameIndex{0};
210  unsigned short m_LastRowNum{0};
211  unsigned short m_LastFrameNum{0};
212  unsigned short m_CurrentRowNum{200};
213  unsigned short m_CurrentFrameNum{2000};
214  unsigned char m_AvailableNodes[14];
215 
216  double m_LocationOffset[14][3];
217 
218  SizeValueType m_NumberOfNodes{0};
219  SizeValueType m_NumberOfCells{0};
220 
221  int m_NodeLimit{2000};
222  int m_CellLimit{4000};
223  int m_ImageWidth{0};
224  int m_ImageHeight{0};
225  int m_ImageDepth{0};
226  int m_ColFlag{0};
227  int m_RowFlag{0};
228  int m_FrameFlag{0};
229  int m_LastRowIndex{0};
230  int m_LastVoxelIndex{0};
231  int m_LastFrameIndex{0};
232 
233  unsigned char m_PointFound{0};
235 
241 };
242 } // end namespace itk
243 
244 #ifndef ITK_MANUAL_INSTANTIATION
245 #include "itkBinaryMask3DMeshSource.hxx"
246 #endif
247 
248 #endif
typename OutputMeshType::MeshTraits OMeshTraits
typename InputImageType::PointType OriginType
typename InputImageType::ConstPointer InputImageConstPointer
Light weight base class for most itk classes.
typename InputImageType::SizeType InputImageSizeType
typename OutputMeshType::PointsContainer PointsContainer
unsigned long SizeValueType
Definition: itkIntTypes.h:83
void SetRegionOfInterest(const RegionType &iRegion)
typename OutputMeshType::PointsContainerPointer PointsContainerPointer
typename OutputMeshType::PointType OPointType
An abstract interface for cells.
typename OutputMeshType::CellsContainerPointer CellsContainerPointer
ImageToMeshFilter is the base class for all process objects that output Mesh data and require image d...
typename InputImageType::IndexType InputImageIndexType
typename OMeshTraits::PixelType OPixelType
SizeValueType IdentifierType
Definition: itkIntTypes.h:87
A multi-dimensional iterator templated over image type that walks a region of pixels.
typename InputImageType::SpacingType SpacingType
typename InputImageType::SizeType SizeType
typename TriCell::SelfAutoPointer TriCellAutoPointer
typename OutputMeshType::CellTraits CellTraits
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename OutputMeshType::Pointer OutputMeshPointer
typename InputImageType::RegionType RegionType
A templated class holding a n-Dimensional covariant vector.
typename InputImageType::Pointer InputImagePointer
typename OutputMeshType::CellsContainer CellsContainer
typename InputImageType::PixelType InputPixelType