ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkGradientRecursiveGaussianImageFilter.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 itkGradientRecursiveGaussianImageFilter_h
19#define itkGradientRecursiveGaussianImageFilter_h
20
23#include "itkImage.h"
24#include "itkCovariantVector.h"
28#include "itkVectorImage.h"
29#include <vector>
30
31namespace itk
32{
54template <
55 typename TInputImage,
56 typename TOutputImage = Image<
58 TInputImage::ImageDimension>>
59class ITK_TEMPLATE_EXPORT GradientRecursiveGaussianImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
60{
61public:
62 ITK_DISALLOW_COPY_AND_MOVE(GradientRecursiveGaussianImageFilter);
63
69
71 using InputImageType = TInputImage;
72 using PixelType = typename TInputImage::PixelType;
75
81
83 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
84
87
92
93
99
101
104
107
110
113
116
118 using OutputImagePointer = typename TOutputImage::Pointer;
119
121 using OutputImageType = TOutputImage;
122 using OutputPixelType = typename OutputImageType::PixelType;
125
127 itkNewMacro(Self);
128
130 itkOverrideGetNameOfClassMacro(GradientRecursiveGaussianImageFilter);
131
133 void
135 void
138
141
144 GetSigma() const;
145
149 void
150 SetNormalizeAcrossScale(bool normalize);
151 itkGetConstMacro(NormalizeAcrossScale, bool);
153
159 void
161
172 itkSetMacro(UseImageDirection, bool);
173 itkGetConstMacro(UseImageDirection, bool);
174 itkBooleanMacro(UseImageDirection);
176
177 // Does not seem to work with wrappings, disabled
178 // itkConceptMacro( InputHasNumericTraitsCheck,
179 // ( Concept::HasNumericTraits< PixelType > ) );
180
181protected:
184 void
185 PrintSelf(std::ostream & os, Indent indent) const override;
186
188 void
189 GenerateData() override;
190
191 // Override since the filter produces the entire dataset
192 void
194
195 void
197
198private:
199 template <typename TValue>
200 void
202 {
203 // To transform Variable length vector we need to convert to and
204 // from the CovariantVectorType
205 const CovariantVectorType gradient(it.Get().GetDataPointer());
206 CovariantVectorType physicalGradient;
207 it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, physicalGradient);
208 it.Set(OutputPixelType(physicalGradient.GetDataPointer(), ImageDimension, false));
209 }
210
211 template <typename T>
212 void
214 {
215 OutputPixelType correctedGradient{};
216 const OutputPixelType & gradient = it.Get();
217
218 const unsigned int nComponents = NumericTraits<OutputPixelType>::GetLength(gradient) / ImageDimension;
219
220 for (unsigned int nc = 0; nc < nComponents; ++nc)
221 {
222 GradientVectorType componentGradient;
223 GradientVectorType correctedComponentGradient;
224 for (unsigned int dim = 0; dim < ImageDimension; ++dim)
225 {
226 componentGradient[dim] =
228 }
229 it.GetImage()->TransformLocalVectorToPhysicalVector(componentGradient, correctedComponentGradient);
230 for (unsigned int dim = 0; dim < ImageDimension; ++dim)
231 {
233 nc * ImageDimension + dim, correctedGradient, correctedComponentGradient[dim]);
234 }
235 }
236 it.Set(correctedGradient);
237 }
238
239 template <template <typename, unsigned int> class P, class T, unsigned int VDimension>
240 void
241 TransformOutputPixel(ImageRegionIterator<Image<P<T, VDimension>, VDimension>> & it)
242 {
243 const OutputPixelType gradient = it.Get();
244 // This uses the more efficient set by reference method
245 it.GetImage()->TransformLocalVectorToPhysicalVector(gradient, it.Value());
246 }
247
248
249 std::vector<GaussianFilterPointer> m_SmoothingFilters{};
252
255
258
261};
262} // end namespace itk
263
264#ifndef ITK_MANUAL_INSTANTIATION
265# include "itkGradientRecursiveGaussianImageFilter.hxx"
266#endif
267
268#endif
A templated class holding a n-Dimensional covariant vector.
Base class for all data objects in ITK.
static ComponentType GetNthComponent(int c, const PixelType &pixel)
static void SetNthComponent(int c, PixelType &pixel, const ComponentType &v)
Simulate a standard C array with copy semantics.
ValueType * GetDataPointer()
typename NumericTraits< PixelType >::RealType RealType
typename NumericTraits< OutputPixelType >::ValueType OutputComponentType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Image< InternalRealType, Self::ImageDimension > RealImageType
typename NumericTraits< InternalRealType >::ValueType InternalScalarRealType
typename NumericTraits< RealType >::FloatType InternalRealType
RecursiveGaussianImageFilter< InputImageType, RealImageType > DerivativeFilterType
void TransformOutputPixel(ImageRegionIterator< VectorImage< TValue, ImageDimension > > &it)
FixedArray< ScalarRealType, Self::ImageDimension > SigmaArrayType
NthElementImageAdaptor< TOutputImage, InternalScalarRealType > OutputImageAdaptorType
typename OutputImageAdaptorType::Pointer OutputImageAdaptorPointer
void SetNormalizeAcrossScale(bool normalize)
void SetSigmaArray(const SigmaArrayType &sigma)
CovariantVector< OutputComponentType, ImageDimension > CovariantVectorType
void TransformOutputPixel(ImageRegionIterator< Image< P< T, VDimension >, VDimension > > &it)
~GradientRecursiveGaussianImageFilter() override=default
typename NumericTraits< PixelType >::ScalarRealType ScalarRealType
CovariantVector< ScalarRealType, ImageDimension > GradientVectorType
RecursiveGaussianImageFilter< RealImageType, RealImageType > GaussianFilterType
void EnlargeOutputRequestedRegion(DataObject *output) override
void SetSigma(ScalarRealType sigma)
void PrintSelf(std::ostream &os, Indent indent) const override
const ImageType * GetImage() const
A multi-dimensional iterator templated over image type that walks a region of pixels.
void Set(const PixelType &value) const
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
Presents an image as being composed of the N-th element of its pixels.
static unsigned int GetLength(const T &)
Base class for computing IIR convolution with an approximation of a Gaussian kernel.
Implements transparent reference counting.
Templated n-dimensional vector image class.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....