ITK  4.1.0
Insight Segmentation and Registration Toolkit
itkDisplacementFieldJacobianDeterminantFilter.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 __itkDisplacementFieldJacobianDeterminantFilter_h
00019 #define __itkDisplacementFieldJacobianDeterminantFilter_h
00020 
00021 #include "itkNeighborhoodIterator.h"
00022 #include "itkImageToImageFilter.h"
00023 #include "itkVector.h"
00024 #include "vnl/vnl_matrix.h"
00025 #include "vnl/vnl_det.h"
00026 
00027 namespace itk
00028 {
00111 template< typename TInputImage,
00112           typename TRealType = float,
00113           typename TOutputImage = Image< TRealType,
00114                                          ::itk::GetImageDimension< TInputImage >::ImageDimension >
00115           >
00116 class ITK_EXPORT DisplacementFieldJacobianDeterminantFilter:
00117   public ImageToImageFilter< TInputImage, TOutputImage >
00118 {
00119 public:
00121   typedef DisplacementFieldJacobianDeterminantFilter      Self;
00122   typedef ImageToImageFilter< TInputImage, TOutputImage > Superclass;
00123   typedef SmartPointer< Self >                            Pointer;
00124   typedef SmartPointer< const Self >                      ConstPointer;
00125 
00127   itkNewMacro(Self);
00128 
00130   itkTypeMacro(DisplacementFieldJacobianDeterminantFilter, ImageToImageFilter);
00131 
00134   typedef typename TOutputImage::PixelType OutputPixelType;
00135   typedef typename TInputImage::PixelType  InputPixelType;
00136 
00138   typedef TInputImage                       InputImageType;
00139   typedef TOutputImage                      OutputImageType;
00140   typedef typename InputImageType::Pointer  InputImagePointer;
00141   typedef typename OutputImageType::Pointer OutputImagePointer;
00142 
00144   itkStaticConstMacro(ImageDimension, unsigned int,
00145                       TOutputImage::ImageDimension);
00146 
00148   itkStaticConstMacro(VectorDimension, unsigned int,
00149                       InputPixelType::Dimension);
00150 
00152   typedef TRealType RealType;
00153   typedef Vector< TRealType, ::itk::GetVectorDimension< InputPixelType >::VectorDimension >
00154   RealVectorType;
00155   typedef Image< RealVectorType, ::itk::GetImageDimension< TInputImage >::ImageDimension >
00156   RealVectorImageType;
00157 
00160   typedef ConstNeighborhoodIterator< RealVectorImageType >   ConstNeighborhoodIteratorType;
00161   typedef typename ConstNeighborhoodIteratorType::RadiusType RadiusType;
00162 
00164   typedef typename Superclass::OutputImageRegionType OutputImageRegionType;
00165 
00174   virtual void GenerateInputRequestedRegion()
00175   throw( InvalidRequestedRegionError );
00176 
00181   void SetUseImageSpacingOn()
00182   { this->SetUseImageSpacing(true); }
00183 
00187   void SetUseImageSpacingOff()
00188   { this->SetUseImageSpacing(false); }
00189 
00192   void SetUseImageSpacing(bool);
00193 
00194   itkGetConstMacro(UseImageSpacing, bool);
00195 
00196   typedef FixedArray< TRealType, ImageDimension > WeightsType;
00197 
00200   void SetDerivativeWeights(const WeightsType &);
00201   itkGetConstReferenceMacro(DerivativeWeights, WeightsType);
00203 
00204 protected:
00205   DisplacementFieldJacobianDeterminantFilter();
00206   virtual ~DisplacementFieldJacobianDeterminantFilter() {}
00207 
00211   void BeforeThreadedGenerateData();
00212 
00225   void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
00226                             ThreadIdType threadId);
00227 
00228   void PrintSelf(std::ostream & os, Indent indent) const;
00229 
00230   typedef typename InputImageType::Superclass ImageBaseType;
00231 
00233   itkGetConstObjectMacro(RealValuedInputImage, ImageBaseType);
00234 
00236   itkGetConstReferenceMacro(NeighborhoodRadius, RadiusType);
00237   itkSetMacro(NeighborhoodRadius, RadiusType);
00239 
00240   virtual TRealType EvaluateAtNeighborhood(const ConstNeighborhoodIteratorType & it) const;
00241 
00243   WeightsType m_DerivativeWeights;
00244 
00247   WeightsType m_HalfDerivativeWeights;
00248 private:
00249   bool m_UseImageSpacing;
00250 
00251   ThreadIdType m_RequestedNumberOfThreads;
00252 
00253   typename ImageBaseType::ConstPointer m_RealValuedInputImage;
00254 
00255   DisplacementFieldJacobianDeterminantFilter(const Self &); //purposely not
00256                                                             // implemented
00257   void operator=(const Self &);                             //purposely not
00258 
00259   // implemented
00260 
00261   RadiusType m_NeighborhoodRadius;
00262 };
00263 } // end namespace itk
00264 
00265 #ifndef ITK_MANUAL_INSTANTIATION
00266 #include "itkDisplacementFieldJacobianDeterminantFilter.hxx"
00267 #endif
00268 
00269 #endif
00270