ITK  4.0.0
Insight Segmentation and Registration Toolkit
itkDiffusionTensor3DReconstructionImageFilter.h
Go to the documentation of this file.
00001 /*=========================================================================
00002  *
00003  *  Copyright Insight Software Consortium
00004  *
00005  *  Licensed under the Apache License, Version 2.0 (the "License");
00006  *  you may not use this file except in compliance with the License.
00007  *  You may obtain a copy of the License at
00008  *
00009  *         http://www.apache.org/licenses/LICENSE-2.0.txt
00010  *
00011  *  Unless required by applicable law or agreed to in writing, software
00012  *  distributed under the License is distributed on an "AS IS" BASIS,
00013  *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
00014  *  See the License for the specific language governing permissions and
00015  *  limitations under the License.
00016  *
00017  *=========================================================================*/
00018 #ifndef __itkDiffusionTensor3DReconstructionImageFilter_h
00019 #define __itkDiffusionTensor3DReconstructionImageFilter_h
00020 
00021 #include "itkImageToImageFilter.h"
00022 #include "itkSpatialObject.h"
00023 #include "itkDiffusionTensor3D.h"
00024 #include "vnl/vnl_matrix.h"
00025 #include "vnl/vnl_vector_fixed.h"
00026 #include "vnl/vnl_matrix_fixed.h"
00027 #include "vnl/algo/vnl_svd.h"
00028 #include "itkVectorContainer.h"
00029 #include "itkVectorImage.h"
00030 
00031 namespace itk
00032 {
00123 template< class TReferenceImagePixelType,
00124           class TGradientImagePixelType = TReferenceImagePixelType,
00125           class TTensorPixelType = double,
00126           class TMaskImageType = Image<unsigned char, 3> >
00127 class ITK_EXPORT DiffusionTensor3DReconstructionImageFilter:
00128   public ImageToImageFilter< Image< TReferenceImagePixelType, 3 >,
00129                              Image< DiffusionTensor3D< TTensorPixelType >, 3 > >
00130 {
00131 public:
00132 
00133   typedef DiffusionTensor3DReconstructionImageFilter Self;
00134   typedef SmartPointer< Self >                       Pointer;
00135   typedef SmartPointer< const Self >                 ConstPointer;
00136   typedef ImageToImageFilter< Image< TReferenceImagePixelType, 3 >,
00137                               Image< DiffusionTensor3D< TTensorPixelType >, 3 > >
00138   Superclass;
00139 
00141   itkNewMacro(Self);
00142 
00144   itkTypeMacro(DiffusionTensor3DReconstructionImageFilter,
00145                ImageToImageFilter);
00146 
00147   typedef TReferenceImagePixelType ReferencePixelType;
00148 
00149   typedef TGradientImagePixelType GradientPixelType;
00150 
00151   typedef DiffusionTensor3D< TTensorPixelType > TensorPixelType;
00152 
00155   typedef typename Superclass::InputImageType ReferenceImageType;
00156 
00157   typedef Image< TensorPixelType, 3 > TensorImageType;
00158 
00159   typedef TensorImageType OutputImageType;
00160 
00161   typedef typename Superclass::OutputImageRegionType
00162   OutputImageRegionType;
00163 
00165   typedef Image< GradientPixelType, 3 > GradientImageType;
00166 
00171   typedef VectorImage< GradientPixelType, 3 > GradientImagesType;
00172 
00174   typedef SpatialObject<3> MaskSpatialObjectType;
00175 
00177   typedef TMaskImageType MaskImageType;
00178 
00180   typedef vnl_matrix_fixed< double, 6, 6 > TensorBasisMatrixType;
00181 
00182   typedef vnl_matrix< double > CoefficientMatrixType;
00183 
00185   typedef vnl_vector_fixed< double, 3 > GradientDirectionType;
00186 
00188   typedef VectorContainer< unsigned int,
00189                            GradientDirectionType >                  GradientDirectionContainerType;
00190 
00192   void AddGradientImage(const GradientDirectionType &, const GradientImageType *image);
00193   const GradientImageType *GetGradientImage(unsigned index) const;
00195 
00202   void SetGradientImage(GradientDirectionContainerType *,
00203                         const GradientImagesType *image);
00204 
00206   void SetReferenceImage(ReferenceImageType *referenceImage)
00207   {
00208     if ( m_GradientImageTypeEnumeration == GradientIsInASingleImage )
00209       {
00210       itkExceptionMacro(<< "Cannot call both methods:"
00211                         << "AddGradientImage and SetGradientImage. Please call only one of them.");
00212       }
00213 
00214     this->ProcessObject::SetNthInput(0, referenceImage);
00215 
00216     m_GradientImageTypeEnumeration = GradientIsInManyImages;
00217   }
00218 
00220   virtual ReferenceImageType * GetReferenceImage()
00221   { return ( static_cast< ReferenceImageType * >( this->ProcessObject::GetInput(0) ) ); }
00222 
00224   virtual GradientDirectionType GetGradientDirection(unsigned int idx) const
00225   {
00226     if ( idx >= m_NumberOfGradientDirections )
00227       {
00228       itkExceptionMacro(<< "Gradient direction " << idx << "does not exist");
00229       }
00230     return m_GradientDirectionContainer->ElementAt(idx + 1);
00231   }
00233 
00235   void SetMaskImage(MaskImageType *maskImage);
00236 
00238   void SetMaskSpatialObject(MaskSpatialObjectType *maskSpatialObject);
00239 
00240 
00244   itkSetMacro(Threshold, ReferencePixelType);
00245   itkGetConstMacro(Threshold, ReferencePixelType);
00247 
00254   itkSetMacro(BValue, TTensorPixelType);
00255 #ifdef GetBValue
00256 #undef GetBValue
00257 #endif
00258   itkGetConstReferenceMacro(BValue, TTensorPixelType);
00260 
00261 #ifdef ITK_USE_CONCEPT_CHECKING
00262 
00263   itkConceptMacro( ReferenceEqualityComparableCheck,
00264                    ( Concept::EqualityComparable< ReferencePixelType > ) );
00265   itkConceptMacro( TensorEqualityComparableCheck,
00266                    ( Concept::EqualityComparable< TensorPixelType > ) );
00267   itkConceptMacro( GradientConvertibleToDoubleCheck,
00268                    ( Concept::Convertible< GradientPixelType, double > ) );
00269   itkConceptMacro( DoubleConvertibleToTensorCheck,
00270                    ( Concept::Convertible< double, TensorPixelType > ) );
00271   itkConceptMacro( GradientReferenceAdditiveOperatorsCheck,
00272                    ( Concept::AdditiveOperators< GradientPixelType, GradientPixelType,
00273                                                  ReferencePixelType > ) );
00274   itkConceptMacro( GradientReferenceAdditiveAndAssignOperatorsCheck,
00275                    ( Concept::AdditiveAndAssignOperators< GradientPixelType,
00276                                                           ReferencePixelType > ) );
00278 
00279   itkConceptMacro( ReferenceOStreamWritableCheck,
00280                    ( Concept::OStreamWritable< ReferencePixelType > ) );
00281   itkConceptMacro( TensorOStreamWritableCheck,
00282                    ( Concept::OStreamWritable< TensorPixelType > ) );
00284 #endif
00285 protected:
00286   DiffusionTensor3DReconstructionImageFilter();
00287   ~DiffusionTensor3DReconstructionImageFilter() {}
00288   void PrintSelf(std::ostream & os, Indent indent) const;
00290 
00291   void ComputeTensorBasis();
00292 
00293   void BeforeThreadedGenerateData();
00294 
00295   void ThreadedGenerateData(const
00296                             OutputImageRegionType & outputRegionForThread, ThreadIdType);
00297 
00300   typedef enum {
00301     GradientIsInASingleImage = 1,
00302     GradientIsInManyImages,
00303     Else
00304     } GradientImageTypeEnumeration;
00305 private:
00306 
00307   /* Tensor basis coeffs */
00308   TensorBasisMatrixType m_TensorBasis;
00309 
00310   CoefficientMatrixType m_BMatrix;
00311 
00313   GradientDirectionContainerType::Pointer m_GradientDirectionContainer;
00314 
00316   unsigned int m_NumberOfGradientDirections;
00317 
00319   unsigned int m_NumberOfBaselineImages;
00320 
00322   ReferencePixelType m_Threshold;
00323 
00325   TTensorPixelType m_BValue;
00326 
00328   GradientImageTypeEnumeration m_GradientImageTypeEnumeration;
00329 
00331   bool m_MaskImagePresent;
00332 };
00333 }
00334 
00335 #ifndef ITK_MANUAL_INSTANTIATION
00336 #include "itkDiffusionTensor3DReconstructionImageFilter.hxx"
00337 #endif
00338 
00339 #endif
00340