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

itkDiffusionTensor3DReconstructionImageFilter.h

Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: itkDiffusionTensor3DReconstructionImageFilter.h,v $
00005   Language:  C++
00006   Date:      $Date: 2005/08/04 13:43:45 $
00007   Version:   $Revision: 1.5 $
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 __itkDiffusionTensor3DReconstructionImageFilter_h_
00018 #define __itkDiffusionTensor3DReconstructionImageFilter_h_
00019 
00020 #include "itkImageToImageFilter.h"
00021 #include "itkDiffusionTensor3D.h"
00022 #include "vnl/vnl_matrix.h"
00023 #include "vnl/vnl_vector_fixed.h"
00024 #include "vnl/vnl_matrix_fixed.h"
00025 #include "vnl/algo/vnl_svd.h"
00026 #include "itkVectorContainer.h"
00027 
00028 namespace itk{
00059 template< class TReferenceImagePixelType, 
00060           class TGradientImagePixelType, class TTensorPixelType=double >
00061 class ITK_EXPORT DiffusionTensor3DReconstructionImageFilter :
00062   public ImageToImageFilter< Image< TReferenceImagePixelType, 3 >, 
00063                              Image< DiffusionTensor3D< TTensorPixelType >, 3 > >
00064 {
00065 
00066 public:
00067 
00068   typedef DiffusionTensor3DReconstructionImageFilter Self;
00069   typedef SmartPointer<Self>                      Pointer;
00070   typedef SmartPointer<const Self>                ConstPointer;
00071   typedef ImageToImageFilter< Image< TReferenceImagePixelType, 3>, 
00072           Image< DiffusionTensor3D< TTensorPixelType >, 3 > >
00073                           Superclass;
00074   
00076   itkNewMacro(Self);  
00077 
00079   itkTypeMacro(DiffusionTensor3DReconstructionImageFilter, 
00080                                                       ImageToImageFilter);
00081  
00082   typedef TReferenceImagePixelType                 ReferencePixelType;
00083 
00084   typedef TGradientImagePixelType                  GradientPixelType;
00085 
00086   typedef DiffusionTensor3D< TTensorPixelType >    TensorPixelType;
00087 
00090   typedef typename Superclass::InputImageType      ReferenceImageType;
00091   
00092   typedef Image< TensorPixelType, 3 >              TensorImageType;
00093   
00094   typedef TensorImageType                          OutputImageType;
00095 
00096   typedef typename Superclass::OutputImageRegionType
00097                                                    OutputImageRegionType;
00098 
00099   typedef Image< GradientPixelType, 3 >            GradientImageType;
00100 
00102   typedef vnl_matrix_fixed< double, 6, 6 >         TensorBasisMatrixType;
00103   
00104   typedef vnl_matrix< double >                     CoefficientMatrixType;
00105 
00107   typedef vnl_vector_fixed< double, 3 >            GradientDirectionType;
00108 
00110   typedef VectorContainer< unsigned int, 
00111           GradientDirectionType >                  GradientDirectionContainerType;
00112   
00113 
00115   void AddGradientImage( const GradientDirectionType &, const GradientImageType *image);
00116 
00118   void SetReferenceImage( ReferenceImageType *referenceImage )
00119     {
00120     this->ProcessObject::SetNthInput( 0, referenceImage );
00121     }
00122     
00124   virtual ReferenceImageType * GetReferenceImage() 
00125   { return ( static_cast< ReferenceImageType *>(this->ProcessObject::GetInput(0)) ); }
00126 
00128   virtual GradientImageType * GetGradientImage( unsigned int idx ) 
00129   {
00130     if( idx >= m_NumberOfGradientDirections )
00131       {
00132       return NULL;
00133       }
00134     return ( static_cast< GradientImageType *>(this->ProcessObject::GetInput(idx+1)) ); 
00135   }
00136 
00138   virtual GradientDirectionType GetGradientDirection( unsigned int idx) const
00139     {
00140     if( idx >= m_NumberOfGradientDirections )
00141       {
00142       itkExceptionMacro( << "Gradient direction " << idx << "does not exist" );
00143       }
00144     return m_GradientDirectionContainer->ElementAt( idx+1 );
00145     }
00146 
00150   itkSetMacro( Threshold, ReferencePixelType );
00151   itkGetMacro( Threshold, ReferencePixelType );
00152 
00153   
00154 protected:
00155   DiffusionTensor3DReconstructionImageFilter();
00156   ~DiffusionTensor3DReconstructionImageFilter() {};
00157   void PrintSelf(std::ostream& os, Indent indent) const;
00158 
00159   void ComputeTensorBasis();
00160   
00161   void BeforeThreadedGenerateData();
00162   void ThreadedGenerateData( const 
00163       OutputImageRegionType &outputRegionForThread, int);
00164   
00165   
00166 private:
00167   
00168   /* Tensor basis coeffs */
00169   TensorBasisMatrixType                             m_TensorBasis;
00170   
00171   CoefficientMatrixType                             m_Coeffs;
00172 
00174   GradientDirectionContainerType::Pointer           m_GradientDirectionContainer;
00175 
00177   unsigned int                                      m_NumberOfGradientDirections;
00178 
00180   ReferencePixelType                                m_Threshold;
00181 };
00182 
00183 }
00184 
00185 #ifndef ITK_MANUAL_INSTANTIATION
00186 #include "itkDiffusionTensor3DReconstructionImageFilter.txx"
00187 #endif
00188 
00189 #endif
00190 

Generated at Wed May 24 23:03:40 2006 for ITK by doxygen 1.3.5 written by Dimitri van Heesch, © 1997-2000