ITK  6.0.0
Insight Toolkit
itkGradientImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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"
25 #include <memory> // For unique_ptr.
26 
27 namespace itk
28 {
29 
30 
31 template <typename TPixelType, unsigned int VImageDimension>
32 class VectorImage;
33 
34 
63 template <typename TInputImage,
64  typename TOperatorValueType = float,
65  typename TOutputValueType = float,
66  typename TOutputImageType =
67  Image<CovariantVector<TOutputValueType, TInputImage::ImageDimension>, TInputImage::ImageDimension>>
68 class ITK_TEMPLATE_EXPORT GradientImageFilter : public ImageToImageFilter<TInputImage, TOutputImageType>
69 {
70 public:
71  ITK_DISALLOW_COPY_AND_MOVE(GradientImageFilter);
75  static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
76  static constexpr unsigned int OutputImageDimension = TOutputImageType::ImageDimension;
77 
79  using InputImageType = TInputImage;
81  using OutputImageType = TOutputImageType;
83 
89 
91  itkNewMacro(Self);
92 
94  itkOverrideGetNameOfClassMacro(GradientImageFilter);
95 
97  using InputPixelType = typename InputImageType::PixelType;
98  using OperatorValueType = TOperatorValueType;
99  using OutputValueType = TOutputValueType;
100  using OutputPixelType = typename OutputImageType::PixelType;
103 
110  void
111  GenerateInputRequestedRegion() override;
112 
118  itkSetMacro(UseImageSpacing, bool);
119  itkGetConstMacro(UseImageSpacing, bool);
120  itkBooleanMacro(UseImageSpacing);
123 #if !defined(ITK_FUTURE_LEGACY_REMOVE)
124 
127  void
128  SetUseImageSpacingOn()
129  {
130  this->SetUseImageSpacing(true);
131  }
132 
136  void
137  SetUseImageSpacingOff()
138  {
139  this->SetUseImageSpacing(false);
140  }
141 #endif
142 
144  void
145  OverrideBoundaryCondition(ImageBoundaryCondition<TInputImage> * boundaryCondition);
146 
147 #ifdef ITK_USE_CONCEPT_CHECKING
148  // Begin concept checking
149  itkConceptMacro(InputConvertibleToOutputCheck, (Concept::Convertible<InputPixelType, OutputValueType>));
150  itkConceptMacro(OutputHasNumericTraitsCheck, (Concept::HasNumericTraits<OutputValueType>));
151  // End concept checking
152 #endif
153 
164  itkSetMacro(UseImageDirection, bool);
165  itkGetConstMacro(UseImageDirection, bool);
166  itkBooleanMacro(UseImageDirection);
169 protected:
170  GradientImageFilter();
171  ~GradientImageFilter() override = default;
172  void
173  PrintSelf(std::ostream & os, Indent indent) const override;
174 
185  void
186  DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
187 
188 
189 private:
190  void
191  GenerateOutputInformation() override;
192 
193  // An overloaded method which may transform the gradient to a
194  // physical vector and converts to the correct output pixel type.
195  template <typename TValue>
196  void
198  {
199  if (this->m_UseImageDirection)
200  {
201  CovariantVectorType physicalGradient;
202  it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, physicalGradient);
203  it.Set(OutputPixelType(physicalGradient.GetDataPointer(), InputImageDimension, false));
204  }
205  else
206  {
207  it.Set(OutputPixelType(gradient.GetDataPointer(), InputImageDimension, false));
208  }
209  }
210 
211  template <typename T>
212  void
214  {
215  // This uses the more efficient set by reference method
216  if (this->m_UseImageDirection)
217  {
218  it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, it.Value());
219  }
220  else
221  {
222  it.Value() = gradient;
223  }
224  }
225 
226 
227  bool m_UseImageSpacing{ true };
228 
229  // flag to take or not the image direction into account
230  // when computing the derivatives.
231  bool m_UseImageDirection{ true };
232 
233  // allow setting the the m_BoundaryCondition
234  std::unique_ptr<ImageBoundaryCondition<TInputImage, TInputImage>> m_BoundaryCondition{
235  std::make_unique<ZeroFluxNeumannBoundaryCondition<TInputImage>>()
236  };
237 };
238 } // end namespace itk
239 
240 #ifndef ITK_MANUAL_INSTANTIATION
241 # include "itkGradientImageFilter.hxx"
242 #endif
243 
244 #endif
itk::GradientImageFilter::InputImageType
TInputImage InputImageType
Definition: itkGradientImageFilter.h:79
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::GradientImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkGradientImageFilter.h:80
itkCovariantVector.h
itk::GradientImageFilter::OutputImageType
TOutputImageType OutputImageType
Definition: itkGradientImageFilter.h:81
itk::VectorImage
Templated n-dimensional vector image class.
Definition: itkImageAlgorithm.h:29
itk::ImageConstIterator::GetImage
const ImageType * GetImage() const
Definition: itkImageConstIterator.h:329
itk::SmartPointer< Self >
itkImageRegionIterator.h
itk::FixedArray::GetDataPointer
ValueType * GetDataPointer()
Definition: itkFixedArray.h:287
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::GradientImageFilter::OutputImageRegionType
typename OutputImageType::RegionType OutputImageRegionType
Definition: itkGradientImageFilter.h:102
itkImageToImageFilter.h
itk::GradientImageFilter
Computes the gradient of an image using directional derivatives.
Definition: itkGradientImageFilter.h:68
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::GradientImageFilter::OutputValueType
TOutputValueType OutputValueType
Definition: itkGradientImageFilter.h:99
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::GradientImageFilter::InputPixelType
typename InputImageType::PixelType InputPixelType
Definition: itkGradientImageFilter.h:97
itk::GradientImageFilter::SetOutputPixel
void SetOutputPixel(ImageRegionIterator< VectorImage< TValue, OutputImageDimension >> &it, CovariantVectorType &gradient)
Definition: itkGradientImageFilter.h:197
itk::GradientImageFilter::SetOutputPixel
void SetOutputPixel(ImageRegionIterator< T > &it, CovariantVectorType &gradient)
Definition: itkGradientImageFilter.h:213
itk::GradientImageFilter::OutputPixelType
typename OutputImageType::PixelType OutputPixelType
Definition: itkGradientImageFilter.h:100
itk::GradientImageFilter::OperatorValueType
TOperatorValueType OperatorValueType
Definition: itkGradientImageFilter.h:98
itk::ImageRegionIterator::Value
PixelType & Value()
Definition: itkImageRegionIterator.h:125
itkZeroFluxNeumannBoundaryCondition.h
itk::GradientImageFilter::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkGradientImageFilter.h:82