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: 2008-10-14 19:20:33 $
00007   Version:   $Revision: 1.13 $
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 #include "itkVectorImage.h"
00028 
00029 namespace itk {
00119 template< class TReferenceImagePixelType, 
00120           class TGradientImagePixelType=TReferenceImagePixelType,
00121           class TTensorPixelType=double >
00122 class ITK_EXPORT DiffusionTensor3DReconstructionImageFilter :
00123   public ImageToImageFilter< Image< TReferenceImagePixelType, 3 >, 
00124                              Image< DiffusionTensor3D< TTensorPixelType >, 3 > >
00125 {
00126 
00127 public:
00128 
00129   typedef DiffusionTensor3DReconstructionImageFilter  Self;
00130   typedef SmartPointer<Self>                          Pointer;
00131   typedef SmartPointer<const Self>                    ConstPointer;
00132   typedef ImageToImageFilter< Image< TReferenceImagePixelType, 3>, 
00133           Image< DiffusionTensor3D< TTensorPixelType >, 3 > >
00134                                                       Superclass;
00135   
00137   itkNewMacro(Self);  
00138 
00140   itkTypeMacro(DiffusionTensor3DReconstructionImageFilter, 
00141                                                       ImageToImageFilter);
00142 
00143   typedef TReferenceImagePixelType                 ReferencePixelType;
00144 
00145   typedef TGradientImagePixelType                  GradientPixelType;
00146 
00147   typedef DiffusionTensor3D< TTensorPixelType >    TensorPixelType;
00148 
00151   typedef typename Superclass::InputImageType      ReferenceImageType;
00152 
00153   typedef Image< TensorPixelType, 3 >              TensorImageType;
00154   
00155   typedef TensorImageType                          OutputImageType;
00156 
00157   typedef typename Superclass::OutputImageRegionType
00158                                                    OutputImageRegionType;
00159 
00161   typedef Image< GradientPixelType, 3 >            GradientImageType;
00162 
00167   typedef VectorImage< GradientPixelType, 3 >      GradientImagesType;
00168 
00170   typedef vnl_matrix_fixed< double, 6, 6 >         TensorBasisMatrixType;
00171 
00172   typedef vnl_matrix< double >                     CoefficientMatrixType;
00173 
00175   typedef vnl_vector_fixed< double, 3 >            GradientDirectionType;
00176 
00178   typedef VectorContainer< unsigned int, 
00179           GradientDirectionType >                  GradientDirectionContainerType;
00180 
00181 
00183   void AddGradientImage( const GradientDirectionType &, const GradientImageType *image);
00184 
00191   void SetGradientImage( GradientDirectionContainerType *, 
00192                                              const GradientImagesType *image);
00193 
00195   void SetReferenceImage( ReferenceImageType *referenceImage )
00196     {
00197     if( m_GradientImageTypeEnumeration == GradientIsInASingleImage)
00198       {
00199       itkExceptionMacro( << "Cannot call both methods:" 
00200       << "AddGradientImage and SetGradientImage. Please call only one of them.");
00201       }
00202 
00203     this->ProcessObject::SetNthInput( 0, referenceImage );
00204 
00205     m_GradientImageTypeEnumeration = GradientIsInManyImages;
00206     }
00207     
00209   virtual ReferenceImageType * GetReferenceImage() 
00210   { return ( static_cast< ReferenceImageType *>(this->ProcessObject::GetInput(0)) ); }
00211 
00213   virtual GradientDirectionType GetGradientDirection( unsigned int idx) const
00214     {
00215     if( idx >= m_NumberOfGradientDirections )
00216       {
00217       itkExceptionMacro( << "Gradient direction " << idx << "does not exist" );
00218       }
00219     return m_GradientDirectionContainer->ElementAt( idx+1 );
00220     }
00222 
00226   itkSetMacro( Threshold, ReferencePixelType );
00227   itkGetMacro( Threshold, ReferencePixelType );
00229 
00230   
00237   itkSetMacro( BValue, TTensorPixelType);
00238 #ifdef GetBValue
00239 #undef GetBValue
00240 #endif
00241   itkGetConstReferenceMacro( BValue, TTensorPixelType);
00243 
00244 #ifdef ITK_USE_CONCEPT_CHECKING
00245 
00246   itkConceptMacro(ReferenceEqualityComparableCheck,
00247     (Concept::EqualityComparable<ReferencePixelType>));
00248   itkConceptMacro(TensorEqualityComparableCheck,
00249     (Concept::EqualityComparable<TensorPixelType>));
00250   itkConceptMacro(GradientConvertibleToDoubleCheck,
00251     (Concept::Convertible<GradientPixelType, double>));
00252   itkConceptMacro(DoubleConvertibleToTensorCheck,
00253     (Concept::Convertible<double, TensorPixelType>));
00254   itkConceptMacro(GradientReferenceAdditiveOperatorsCheck,
00255     (Concept::AdditiveOperators<GradientPixelType, GradientPixelType,
00256                                 ReferencePixelType>));
00257   itkConceptMacro(ReferenceOStreamWritableCheck,
00258     (Concept::OStreamWritable<ReferencePixelType>));
00259   itkConceptMacro(TensorOStreamWritableCheck,
00260     (Concept::OStreamWritable<TensorPixelType>));
00261 
00263 #endif
00264 
00265 protected:
00266   DiffusionTensor3DReconstructionImageFilter();
00267   ~DiffusionTensor3DReconstructionImageFilter() {};
00268   void PrintSelf(std::ostream& os, Indent indent) const;
00269 
00270   void ComputeTensorBasis();
00271   
00272   void BeforeThreadedGenerateData();
00273   void ThreadedGenerateData( const 
00274       OutputImageRegionType &outputRegionForThread, int);
00275   
00278   typedef enum
00279     {
00280     GradientIsInASingleImage = 1,
00281     GradientIsInManyImages,
00282     Else
00283     } GradientImageTypeEnumeration;
00284 
00285 private:
00286   
00287   /* Tensor basis coeffs */
00288   TensorBasisMatrixType                             m_TensorBasis;
00289   
00290   CoefficientMatrixType                             m_BMatrix;
00291 
00293   GradientDirectionContainerType::Pointer           m_GradientDirectionContainer;
00294 
00296   unsigned int                                      m_NumberOfGradientDirections;
00297 
00299   unsigned int                                      m_NumberOfBaselineImages;
00300 
00302   ReferencePixelType                                m_Threshold;
00303 
00305   TTensorPixelType                                  m_BValue;
00306 
00308   GradientImageTypeEnumeration                      m_GradientImageTypeEnumeration;
00309 };
00310 
00311 }
00312 
00313 #ifndef ITK_MANUAL_INSTANTIATION
00314 #include "itkDiffusionTensor3DReconstructionImageFilter.txx"
00315 #endif
00316 
00317 #endif
00318 

Generated at Sat Feb 28 12:19:36 2009 for ITK by doxygen 1.5.6 written by Dimitri van Heesch, © 1997-2000