ITK  4.9.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 BinaryMask3DMeshSource:public ImageToMeshFilter< TInputImage, TOutputMesh >
72 {
73 public:
79 
81  itkNewMacro(Self);
82 
85 
87  typedef TOutputMesh OutputMeshType;
88  typedef typename OutputMeshType::MeshTraits OMeshTraits;
89  typedef typename OutputMeshType::PointType OPointType;
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;
114  typedef typename InputImageType::PointType OriginType;
115  typedef typename InputImageType::RegionType RegionType;
116  typedef typename InputImageType::SizeType SizeType;
117 
119  typedef typename InputImageType::IndexType InputImageIndexType;
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;
141  this->Modified();
142  }
143  }
144 
145  itkGetConstReferenceMacro(RegionOfInterest, RegionType);
146 
147 protected:
150  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
151 
152  void GenerateData() ITK_OVERRIDE;
153 
154 
157 
158  virtual void GenerateOutputInformation() ITK_OVERRIDE {} // do nothing ITK_OVERRIDE
159 
160 private:
161  BinaryMask3DMeshSource(const Self &) ITK_DELETE_FUNCTION;
162  void operator=(const Self &) ITK_DELETE_FUNCTION;
163 
165 
166  void CreateMesh();
167 
168  void XFlip(unsigned char *tp); // 7 kinds of transformation
169 
170  void YFlip(unsigned char *tp);
171 
172  void ZFlip(unsigned char *tp);
173 
174  void XRotation(unsigned char *tp);
175 
176  void YRotation(unsigned char *tp);
177 
178  void ZRotation(unsigned char *tp);
179 
180  void inverse(unsigned char *tp);
181 
182  void InitializeLUT(); // initialize the look up table before the mesh
183  // construction
184 
185  void AddCells(unsigned char celltype, unsigned char celltran, int index);
186 
187  void AddNodes(int index,
188  unsigned char *nodesid,
189  IdentifierType *globalnodesid,
190  IdentifierType **currentrowtmp,
191  IdentifierType **currentframetmp);
192 
193  void CellTransfer(unsigned char *nodesid, unsigned char celltran);
194 
195  IdentifierType SearchThroughLastRow(int index, int start, int end);
196 
197  IdentifierType SearchThroughLastFrame(int index, int start, int end);
198 
199  unsigned char m_LUT[256][2]; // the two lookup tables
200 
203 
208 
209  unsigned short m_CurrentRowIndex;
210  unsigned short m_CurrentFrameIndex;
211  unsigned short m_LastRowNum;
212  unsigned short m_LastFrameNum;
213  unsigned short m_CurrentRowNum;
214  unsigned short m_CurrentFrameNum;
215  unsigned char m_AvailableNodes[14];
216 
217  double m_LocationOffset[14][3];
218 
221 
233 
234  unsigned char m_PointFound;
236 
242 };
243 } // end namespace itk
244 
245 #ifndef ITK_MANUAL_INSTANTIATION
246 #include "itkBinaryMask3DMeshSource.hxx"
247 #endif
248 
249 #endif
OutputMeshType::MeshTraits OMeshTraits
void YRotation(unsigned char *tp)
void PrintSelf(std::ostream &os, Indent indent) const override
Light weight base class for most itk classes.
OutputMeshType::CellTraits CellTraits
IdentifierType SearchThroughLastRow(int index, int start, int end)
ImageToMeshFilter< TInputImage, TOutputMesh > Superclass
SmartPointer< const Self > ConstPointer
void AddNodes(int index, unsigned char *nodesid, IdentifierType *globalnodesid, IdentifierType **currentrowtmp, IdentifierType **currentframetmp)
void XFlip(unsigned char *tp)
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
void CellTransfer(unsigned char *nodesid, unsigned char celltran)
CovariantVector< double, 2 > doubleVector
virtual void GenerateOutputInformation() override
TriangleCell< TCellInterface > TriCell
InputImageType::RegionType RegionType
void XRotation(unsigned char *tp)
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
virtual void SetInput(const InputImageType *inputImage)
OutputMeshType::Pointer OutputMeshPointer
void ZRotation(unsigned char *tp)
void AddCells(unsigned char celltype, unsigned char celltran, int index)
A multi-dimensional iterator templated over image type that walks a region of pixels.
CovariantVector< int, 2 > intVector
OutputMeshType::PointType OPointType
void inverse(unsigned char *tp)
virtual void Modified() const
OutputMeshType::CellsContainer CellsContainer
void ZFlip(unsigned char *tp)
void GenerateData() override
InputImageType::Pointer InputImagePointer
InputImageType::SpacingType SpacingType
InputImageType::PixelType InputPixelType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
void SetInput(unsigned int idx, const InputImageType *input)
void YFlip(unsigned char *tp)
InputImageType::PointType OriginType
InputImageType::SizeType InputImageSizeType
InputImageType::ConstPointer InputImageConstPointer
IdentifierType SearchThroughLastFrame(int index, int start, int end)
A templated class holding a n-Dimensional covariant vector.
OutputMeshType::PointsContainerPointer PointsContainerPointer
OutputMeshType::CellsContainerPointer CellsContainerPointer