Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

itkTriangleMeshToBinaryImageFilter.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkTriangleMeshToBinaryImageFilter.h,v $
00005   Language:  C++
00006   Date:      $Date: 2010-06-09 08:46:19 $
00007   Version:   $Revision: 1.12 $
00008 
00009   Copyright (c) Insight Software Consortium. All rights reserved.
00010   See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
00011 
00012      This software is distributed WITHOUT ANY WARRANTY; without even 
00013      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
00014      PURPOSE.  See the above copyright notices for more information.
00015 
00016 =========================================================================*/
00017 #ifndef __itkTriangleMeshToBinaryImageFilter_h
00018 #define __itkTriangleMeshToBinaryImageFilter_h
00019 
00020 #include "itkImageSource.h"
00021 
00022 #include <itkMesh.h>
00023 #include <itkLineCell.h>
00024 #include <itkPolygonCell.h>
00025 #include <itkVertexCell.h>
00026 #include <itkMapContainer.h>
00027 #include <itkVectorContainer.h>
00028 #include <itkVector.h>
00029 #include <itkPoint.h>
00030 #include <itkArray.h>
00031 #include "itkAutomaticTopologyMeshSource.h"
00032 #include "itkPointSet.h"
00033 
00034 #include <vector>
00035 
00036 namespace itk
00037 {
00038 class Point1D  
00039 {
00040 public:
00041   double m_X;
00042   int    m_Sign;
00043      
00044   Point1D(){};
00045   Point1D(const double p, const int s)
00046     {
00047     m_X = p;
00048     m_Sign = s;
00049     }
00050   Point1D(const Point1D &point)
00051     {
00052     m_X = point.m_X;
00053     m_Sign = point.m_Sign;
00054     }
00055   double getX() const
00056     {
00057     return m_X;
00058     }
00059   int  getSign() const
00060     {
00061     return m_Sign;
00062     }
00063   
00064 };
00065 
00072 template <class TInputMesh, class TOutputImage>
00073 class ITK_EXPORT TriangleMeshToBinaryImageFilter : public ImageSource<TOutputImage>
00074 {
00075 public:
00077   typedef TriangleMeshToBinaryImageFilter  Self;
00078   typedef ImageSource<TOutputImage>        Superclass;
00079   typedef SmartPointer<Self>               Pointer;
00080   typedef SmartPointer<const Self>         ConstPointer;
00081 
00082   typedef typename TOutputImage::IndexType        IndexType;
00083   typedef typename TOutputImage::SizeType         SizeType;
00084   typedef TOutputImage                            OutputImageType;
00085   typedef typename OutputImageType::Pointer       OutputImagePointer;
00086   typedef typename OutputImageType::ValueType     ValueType;
00087   typedef typename OutputImageType::SpacingType   SpacingType;
00088   typedef typename OutputImageType::DirectionType DirectionType;
00089   
00091   itkNewMacro(Self);
00092 
00094   itkTypeMacro(TriangleMeshToBinaryImageFilter,ImageSource);
00095 
00097   typedef typename Superclass::OutputImageRegionType OutputImageRegionType;
00098 
00100   typedef TInputMesh                                     InputMeshType;
00101   typedef typename InputMeshType::Pointer                InputMeshPointer;
00102   typedef typename InputMeshType::PointType              InputPointType;
00103   typedef typename InputMeshType::PixelType              InputPixelType;
00104   typedef typename InputMeshType::MeshTraits::CellTraits InputCellTraitsType;
00105   typedef typename InputMeshType::CellType               CellType;
00106   typedef typename InputMeshType::CellsContainerPointer  CellsContainerPointer;
00107   typedef typename InputMeshType::CellsContainerIterator CellsContainerIterator;
00108 
00109   typedef typename InputMeshType::PointsContainer InputPointsContainer;
00110   typedef typename InputPointsContainer::Pointer  InputPointsContainerPointer;
00111   typedef typename InputPointsContainer::Iterator InputPointsContainerIterator;
00112   
00113   typedef itk::PointSet < double , 3 >           PointSetType;
00114   typedef typename PointSetType::PointsContainer PointsContainer;
00115   
00116 
00117   typedef itk::Point < double, 3> PointType;
00118   typedef itk::Point < double, 2> Point2DType;
00119   
00120   typedef itk::Array<double> DoubleArrayType;
00121 
00122   typedef std::vector< Point1D >               Point1DVector;
00123   typedef std::vector< std::vector< Point1D> > Point1DArray;
00124 
00125   typedef std::vector< Point2DType >                Point2DVector;
00126   typedef std::vector< std::vector< Point2DType > > Point2DArray;
00127 
00128   typedef std::vector< PointType >                PointVector;
00129   typedef std::vector< std::vector< PointType > > PointArray;
00130 
00131   typedef std::vector< int > StencilIndexVector;
00136   itkSetMacro(Spacing, SpacingType);
00137   virtual void SetSpacing( const double spacing[3] );
00138   virtual void SetSpacing( const float spacing[3] );
00139   itkGetConstReferenceMacro(Spacing, SpacingType);
00141 
00145   itkSetMacro(Direction, DirectionType);
00146   itkGetConstMacro(Direction, DirectionType);
00148 
00154   itkSetMacro(InsideValue, ValueType);
00155   itkGetConstMacro(InsideValue, ValueType);
00157 
00163   itkSetMacro(OutsideValue, ValueType);
00164   itkGetConstMacro(OutsideValue, ValueType);
00166 
00171   itkSetMacro(Origin, PointType);
00172   virtual void SetOrigin( const double origin[3] );
00173   virtual void SetOrigin( const float origin[3] );
00174   itkGetConstReferenceMacro(Origin, PointType);
00176 
00178   itkSetMacro(Index,IndexType);
00179   itkGetConstMacro(Index,IndexType);
00181 
00183   itkSetMacro(Size,SizeType);
00184   itkGetConstMacro(Size,SizeType);
00186 
00188   void SetInput(InputMeshType *input);
00189 
00190   void SetInfoImage(OutputImageType *InfoImage)
00191     {
00192     if (InfoImage != m_InfoImage)
00193       {
00194       this->Modified();
00195       m_InfoImage = InfoImage;
00196       }
00197     }
00198 
00200   InputMeshType * GetInput(void);
00201   InputMeshType * GetInput(unsigned int idx);
00203 
00204   /* Set the tolerance for doing spatial searches of the polydata. */
00205   itkSetMacro(Tolerance, double);
00206   itkGetConstMacro(Tolerance, double);
00207     
00208 protected:
00209   TriangleMeshToBinaryImageFilter();
00210   ~TriangleMeshToBinaryImageFilter();
00211 
00212   virtual void GenerateOutputInformation(){}; // do nothing
00213   virtual void GenerateData();
00214 
00215   virtual void RasterizeTriangles();
00216   static int PolygonToImageRaster( PointVector coords, Point1DArray & zymatrix, int extent[6]);
00217   
00218   OutputImageType *m_InfoImage;
00219   IndexType        m_Index;
00220   SizeType         m_Size;
00221   SpacingType      m_Spacing;
00222   PointType        m_Origin; //start value
00223   double           m_Tolerance;
00224   ValueType        m_InsideValue;
00225   ValueType        m_OutsideValue;
00226   DirectionType    m_Direction;
00227 
00228   StencilIndexVector m_StencilIndex;
00229 
00230   virtual void PrintSelf(std::ostream& os, Indent indent) const;
00231 
00232 private:
00233   TriangleMeshToBinaryImageFilter(const Self&); //purposely not implemented
00234   void operator=(const Self&); //purposely not implemented
00235 
00236   static bool ComparePoints2D(Point2DType a, Point2DType b);
00237   static bool ComparePoints1D(Point1D a, Point1D b);
00238 
00239 };
00240 
00241 } // end namespace itk
00242 
00243 #ifndef ITK_MANUAL_INSTANTIATION
00244 #include "itkTriangleMeshToBinaryImageFilter.txx"
00245 #endif
00246 
00247 #endif
00248 

Generated at Mon Jul 12 2010 20:03:55 for ITK by doxygen 1.7.1 written by Dimitri van Heesch, © 1997-2000