ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkGradientMagnitudeImageFilter.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 itkGradientMagnitudeImageFilter_h
19#define itkGradientMagnitudeImageFilter_h
20
22
23namespace itk
24{
41template <typename TInputImage, typename TOutputImage>
42class ITK_TEMPLATE_EXPORT GradientMagnitudeImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
43{
44public:
45 ITK_DISALLOW_COPY_AND_MOVE(GradientMagnitudeImageFilter);
46
52
54 itkNewMacro(Self);
55
57 itkOverrideGetNameOfClassMacro(GradientMagnitudeImageFilter);
58
61 using OutputPixelType = typename TOutputImage::PixelType;
62 using InputPixelType = typename TInputImage::PixelType;
64
67 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
68
70 using InputImageType = TInputImage;
71 using OutputImageType = TOutputImage;
72 using InputImagePointer = typename InputImageType::Pointer;
73 using OutputImagePointer = typename OutputImageType::Pointer;
74
77
86 void
88
93 itkSetMacro(UseImageSpacing, bool);
94 itkGetConstMacro(UseImageSpacing, bool);
95 itkBooleanMacro(UseImageSpacing);
97
98#if !defined(ITK_FUTURE_LEGACY_REMOVE)
102 void
103 SetUseImageSpacingOn()
104 {
105 this->SetUseImageSpacing(true);
106 }
107
111 void
112 SetUseImageSpacingOff()
113 {
114 this->SetUseImageSpacing(false);
115 }
116#endif
117
118 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<InputPixelType>));
119
120protected:
122
123 ~GradientMagnitudeImageFilter() override = default;
124
136 void
137 DynamicThreadedGenerateData(const OutputImageRegionType & outputRegionForThread) override;
138
139
140 void
141 PrintSelf(std::ostream &, Indent) const override;
142
143private:
144 bool m_UseImageSpacing{ true };
145};
146} // end namespace itk
147
148#ifndef ITK_MANUAL_INSTANTIATION
149# include "itkGradientMagnitudeImageFilter.hxx"
150#endif
151
152#endif
typename NumericTraits< InputPixelType >::RealType RealType
typename OutputImageType::Pointer OutputImagePointer
typename TOutputImage::PixelType OutputPixelType
~GradientMagnitudeImageFilter() override=default
typename InputImageType::Pointer InputImagePointer
ImageToImageFilter< TInputImage, TOutputImage > Superclass
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
void GenerateInputRequestedRegion() override
void PrintSelf(std::ostream &, Indent) const override
typename TInputImage::PixelType InputPixelType
typename OutputImageType::RegionType OutputImageRegionType
Control indentation during Print() invocation.
Definition itkIndent.h:50
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....