ITK  4.13.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:
79 
81  itkNewMacro(Self);
82 
85 
87  typedef TOutputMesh OutputMeshType;
88  typedef typename OutputMeshType::MeshTraits OMeshTraits;
90  typedef typename OMeshTraits::PixelType OPixelType;
91 
93  typedef typename OutputMeshType::Pointer OutputMeshPointer;
94  typedef typename OutputMeshType::CellTraits CellTraits;
95  typedef typename OutputMeshType::PointsContainerPointer PointsContainerPointer;
96  typedef typename OutputMeshType::PointsContainer PointsContainer;
97  typedef typename OutputMeshType::CellsContainerPointer CellsContainerPointer;
98  typedef typename OutputMeshType::CellsContainer CellsContainer;
101 
106  typedef typename TriCell::SelfAutoPointer TriCellAutoPointer;
107 
109  typedef TInputImage InputImageType;
110  typedef typename InputImageType::Pointer InputImagePointer;
111  typedef typename InputImageType::ConstPointer InputImageConstPointer;
112  typedef typename InputImageType::PixelType InputPixelType;
113  typedef typename InputImageType::SpacingType SpacingType;
115  typedef typename InputImageType::RegionType RegionType;
117 
120 
122 
125 
126  itkSetMacro(ObjectValue, InputPixelType);
127 
128  itkGetConstMacro(NumberOfNodes, SizeValueType);
129  itkGetConstMacro(NumberOfCells, SizeValueType);
130 
132  using Superclass::SetInput;
133  virtual void SetInput(const InputImageType *inputImage);
134 
135  void SetRegionOfInterest( const RegionType & iRegion )
136  {
137  if( iRegion != m_RegionOfInterest )
138  {
139  this->m_RegionOfInterest = iRegion;
140  this->m_RegionOfInterestProvidedByUser = true;
141  this->Modified();
142  }
143  }
144 
145  itkGetConstReferenceMacro(RegionOfInterest, RegionType);
146 
147 protected:
149  ~BinaryMask3DMeshSource() ITK_OVERRIDE;
150  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
151 
152  void GenerateData() ITK_OVERRIDE;
153 
154 
155  bool m_RegionOfInterestProvidedByUser;
156  RegionType m_RegionOfInterest;
157 
158  virtual void GenerateOutputInformation() ITK_OVERRIDE {} // do nothing ITK_OVERRIDE
159 
160 private:
161  ITK_DISALLOW_COPY_AND_ASSIGN(BinaryMask3DMeshSource);
162 
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 
207 
208  unsigned short m_CurrentRowIndex;
209  unsigned short m_CurrentFrameIndex;
210  unsigned short m_LastRowNum;
211  unsigned short m_LastFrameNum;
212  unsigned short m_CurrentRowNum;
213  unsigned short m_CurrentFrameNum;
214  unsigned char m_AvailableNodes[14];
215 
216  double m_LocationOffset[14][3];
217 
220 
232 
233  unsigned char m_PointFound;
235 
241 };
242 } // end namespace itk
243 
244 #ifndef ITK_MANUAL_INSTANTIATION
245 #include "itkBinaryMask3DMeshSource.hxx"
246 #endif
247 
248 #endif
OutputMeshType::MeshTraits OMeshTraits
Light weight base class for most itk classes.
OutputMeshType::CellTraits CellTraits
ImageToMeshFilter< TInputImage, TOutputMesh > Superclass
SmartPointer< const Self > ConstPointer
void SetRegionOfInterest(const RegionType &iRegion)
ImageRegionConstIterator< InputImageType > InputImageIterator
InputImageType::IndexType InputImageIndexType
TriCell::SelfAutoPointer TriCellAutoPointer
CellInterface< OPixelType, CellTraits > TCellInterface
An abstract interface for cells.
unsigned long SizeValueType
Definition: itkIntTypes.h:143
CovariantVector< double, 2 > doubleVector
TriangleCell< TCellInterface > TriCell
InputImageType::RegionType RegionType
ImageToMeshFilter is the base class for all process objects that output Mesh data and require image d...
OutputMeshType::PointsContainer PointsContainer
SizeValueType IdentifierType
Definition: itkIntTypes.h:147
OutputMeshType::Pointer OutputMeshPointer
A multi-dimensional iterator templated over image type that walks a region of pixels.
CovariantVector< int, 2 > intVector
OutputMeshType::PointType OPointType
OutputMeshType::CellsContainer CellsContainer
InputImageType::Pointer InputImagePointer
InputImageType::SpacingType SpacingType
InputImageType::PixelType InputPixelType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
InputImageType::PointType OriginType
InputImageType::SizeType InputImageSizeType
InputImageType::ConstPointer InputImageConstPointer
A templated class holding a n-Dimensional covariant vector.
OutputMeshType::PointsContainerPointer PointsContainerPointer
OutputMeshType::CellsContainerPointer CellsContainerPointer