ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkGradientImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkGradientImageFilter_h
19 #define itkGradientImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkCovariantVector.h"
23 #include "itkImageRegionIterator.h"
24 
25 namespace itk
26 {
27 
28 
29 template <typename TPixelType, unsigned int VImageDimension > class VectorImage;
30 
31 
54 template< typename TInputImage,
55  typename TOperatorValueType = float,
56  typename TOutputValueType = float,
57  typename TOutputImageType = Image< CovariantVector< TOutputValueType,
58  TInputImage::ImageDimension >,
59  TInputImage::ImageDimension > >
61  public ImageToImageFilter< TInputImage, TOutputImageType >
62 {
63 public:
65  itkStaticConstMacro(InputImageDimension, unsigned int,
66  TInputImage::ImageDimension);
67  itkStaticConstMacro(OutputImageDimension, unsigned int,
68  TInputImage::ImageDimension);
70 
73 
75  typedef TInputImage InputImageType;
76  typedef typename InputImageType::Pointer InputImagePointer;
77  typedef TOutputImageType OutputImageType;
78  typedef typename OutputImageType::Pointer OutputImagePointer;
79 
84 
86  itkNewMacro(Self);
87 
90 
92  typedef typename InputImageType::PixelType InputPixelType;
93  typedef TOperatorValueType OperatorValueType;
94  typedef TOutputValueType OutputValueType;
95  typedef typename OutputImageType::PixelType OutputPixelType;
96  typedef CovariantVector<
97  OutputValueType, itkGetStaticConstMacro(OutputImageDimension) >
99  typedef typename OutputImageType::RegionType OutputImageRegionType;
100 
107  virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
108 
112  { this->SetUseImageSpacing(true); }
113 
117  { this->SetUseImageSpacing(false); }
118 
121  itkSetMacro(UseImageSpacing, bool);
122  itkGetConstMacro(UseImageSpacing, bool);
124 
125 #ifdef ITK_USE_CONCEPT_CHECKING
126  // Begin concept checking
127  itkConceptMacro( InputConvertibleToOutputCheck,
129  itkConceptMacro( OutputHasNumericTraitsCheck,
131  // End concept checking
132 #endif
133 
144  itkSetMacro(UseImageDirection, bool);
145  itkGetConstMacro(UseImageDirection, bool);
146  itkBooleanMacro(UseImageDirection);
148 
149 protected:
151  virtual ~GradientImageFilter();
152  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
153 
164  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
165  ThreadIdType threadId) ITK_OVERRIDE;
166 
167 private:
168  GradientImageFilter(const Self &); //purposely not implemented
169  void operator=(const Self &); //purposely not implemented
170 
171  virtual void GenerateOutputInformation() ITK_OVERRIDE;
172 
173  // An overloaded method which may transform the gradient to a
174  // physical vector and converts to the correct output pixel type.
175  template <typename TValue>
177  {
178  if ( this->m_UseImageDirection )
179  {
180  CovariantVectorType physicalGradient;
181  it.GetImage()->TransformLocalVectorToPhysicalVector( gradient, physicalGradient );
182  it.Set( OutputPixelType( physicalGradient.GetDataPointer(), InputImageDimension, false ) );
183  }
184  else
185  {
186  it.Set( OutputPixelType( gradient.GetDataPointer(), InputImageDimension, false ) );
187  }
188  }
189 
190  template <typename T >
192  {
193  // This uses the more efficient set by reference method
194  if ( this->m_UseImageDirection )
195  {
196  it.GetImage()->TransformLocalVectorToPhysicalVector( gradient, it.Value() );
197  }
198  else
199  {
200  it.Value() = gradient;
201  }
202  }
203 
204 
206 
207  // flag to take or not the image direction into account
208  // when computing the derivatives.
210 };
211 } // end namespace itk
212 
213 #ifndef ITK_MANUAL_INSTANTIATION
214 #include "itkGradientImageFilter.hxx"
215 #endif
216 
217 #endif
Light weight base class for most itk classes.
InputImageType::Pointer InputImagePointer
TOperatorValueType OperatorValueType
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId) override
InputImageType::PixelType InputPixelType
Templated n-dimensional vector image class.
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void GenerateOutputInformation() override
OutputImageType::RegionType OutputImageRegionType
void SetOutputPixel(ImageRegionIterator< VectorImage< TValue, OutputImageDimension > > &it, CovariantVectorType &gradient)
OutputImageType::Pointer OutputImagePointer
SmartPointer< const Self > ConstPointer
void SetOutputPixel(ImageRegionIterator< T > &it, CovariantVectorType &gradient)
CovariantVector< OutputValueType, itkGetStaticConstMacro(OutputImageDimension) > CovariantVectorType
void operator=(const Self &)
static const unsigned int OutputImageDimension
SmartPointer< Self > Pointer
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
virtual void SetUseImageSpacing(bool _arg)
OutputImageType::PixelType OutputPixelType
ImageToImageFilter< InputImageType, OutputImageType > Superclass
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
virtual void GenerateInputRequestedRegion() override
ValueType * GetDataPointer()
#define itkConceptMacro(name, concept)
static const unsigned int InputImageDimension
A templated class holding a n-Dimensional covariant vector.
A multi-dimensional iterator templated over image type that walks a region of pixels.
const ImageType * GetImage() const
Computes the gradient of an image using directional derivatives.