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: 2009-04-25 12:28:11 $
00007   Version:   $Revision: 1.11 $
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   
00090   itkNewMacro(Self);
00091 
00093   itkTypeMacro(TriangleMeshToBinaryImageFilter,ImageSource);
00094 
00096   typedef typename Superclass::OutputImageRegionType OutputImageRegionType;
00097 
00099   typedef TInputMesh                                     InputMeshType;
00100   typedef typename InputMeshType::Pointer                InputMeshPointer;
00101   typedef typename InputMeshType::PointType              InputPointType;
00102   typedef typename InputMeshType::PixelType              InputPixelType;
00103   typedef typename InputMeshType::MeshTraits::CellTraits InputCellTraitsType;
00104   typedef typename InputMeshType::CellType               CellType;
00105   typedef typename InputMeshType::CellsContainerPointer  CellsContainerPointer;
00106   typedef typename InputMeshType::CellsContainerIterator CellsContainerIterator;
00107 
00108   typedef typename InputMeshType::PointsContainer InputPointsContainer;
00109   typedef typename InputPointsContainer::Pointer  InputPointsContainerPointer;
00110   typedef typename InputPointsContainer::Iterator InputPointsContainerIterator;
00111   
00112   typedef itk::PointSet < double , 3 >           PointSetType;
00113   typedef typename PointSetType::PointsContainer PointsContainer;
00114   
00115 
00116   typedef itk::Point < double, 3> PointType;
00117   typedef itk::Point < double, 2> Point2DType;
00118   
00119   typedef itk::Array<double> DoubleArrayType;
00120 
00121   typedef std::vector< Point1D >               Point1DVector;
00122   typedef std::vector< std::vector< Point1D> > Point1DArray;
00123 
00124   typedef std::vector< Point2DType >                Point2DVector;
00125   typedef std::vector< std::vector< Point2DType > > Point2DArray;
00126 
00127   typedef std::vector< PointType >                PointVector;
00128   typedef std::vector< std::vector< PointType > > PointArray;
00129 
00130   typedef std::vector< int > StencilIndexVector;
00135   itkSetMacro(Spacing, SpacingType);
00136   virtual void SetSpacing( const double spacing[3] );
00137   virtual void SetSpacing( const float spacing[3] );
00138   itkGetConstReferenceMacro(Spacing, SpacingType);
00140 
00146   itkSetMacro(InsideValue, ValueType);
00147   itkGetConstMacro(InsideValue, ValueType);
00149 
00155   itkSetMacro(OutsideValue, ValueType);
00156   itkGetConstMacro(OutsideValue, ValueType);
00158 
00163   itkSetMacro(Origin, PointType);
00164   virtual void SetOrigin( const double origin[3] );
00165   virtual void SetOrigin( const float origin[3] );
00166   itkGetConstReferenceMacro(Origin, PointType);
00168 
00170   itkSetMacro(Index,IndexType);
00171   itkGetConstMacro(Index,IndexType);
00173 
00175   itkSetMacro(Size,SizeType);
00176   itkGetConstMacro(Size,SizeType);
00178 
00180   void SetInput(InputMeshType *input);
00181 
00182   void SetInfoImage(OutputImageType *InfoImage)
00183     {
00184     if (InfoImage != m_InfoImage)
00185       {
00186       this->Modified();
00187       m_InfoImage = InfoImage;
00188       }
00189     }
00190 
00192   InputMeshType * GetInput(void);
00193   InputMeshType * GetInput(unsigned int idx);
00195 
00196   /* Set the tolerance for doing spatial searches of the polydata. */
00197   itkSetMacro(Tolerance, double);
00198   itkGetConstMacro(Tolerance, double);
00199     
00200 protected:
00201   TriangleMeshToBinaryImageFilter();
00202   ~TriangleMeshToBinaryImageFilter();
00203 
00204   virtual void GenerateOutputInformation(){}; // do nothing
00205   virtual void GenerateData();
00206 
00207   virtual void RasterizeTriangles();
00208   static int PolygonToImageRaster( PointVector coords, Point1DArray & zymatrix, int extent[6]);
00209   
00210   OutputImageType *m_InfoImage;
00211   IndexType        m_Index;
00212   SizeType         m_Size;
00213   SpacingType      m_Spacing;
00214   PointType        m_Origin; //start value
00215   double           m_Tolerance;
00216   ValueType        m_InsideValue;
00217   ValueType        m_OutsideValue;
00218 
00219   StencilIndexVector m_StencilIndex;
00220 
00221   virtual void PrintSelf(std::ostream& os, Indent indent) const;
00222 
00223 private:
00224   TriangleMeshToBinaryImageFilter(const Self&); //purposely not implemented
00225   void operator=(const Self&); //purposely not implemented
00226 
00227   static bool ComparePoints2D(Point2DType a, Point2DType b);
00228   static bool ComparePoints1D(Point1D a, Point1D b);
00229 
00230 };
00231 
00232 } // end namespace itk
00233 
00234 #ifndef ITK_MANUAL_INSTANTIATION
00235 #include "itkTriangleMeshToBinaryImageFilter.txx"
00236 #endif
00237 
00238 #endif
00239 

Generated at Fri Apr 16 19:50:28 2010 for ITK by doxygen 1.6.1 written by Dimitri van Heesch, © 1997-2000