ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkDiscreteGaussianImageFilter.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 itkDiscreteGaussianImageFilter_h
19#define itkDiscreteGaussianImageFilter_h
20
21#include "itkGaussianOperator.h"
23#include "itkImage.h"
25
26namespace itk
27{
60
61template <typename TInputImage, typename TOutputImage = TInputImage>
62class ITK_TEMPLATE_EXPORT DiscreteGaussianImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
63{
64public:
65 ITK_DISALLOW_COPY_AND_MOVE(DiscreteGaussianImageFilter);
66
72
74 itkNewMacro(Self);
75
77 itkOverrideGetNameOfClassMacro(DiscreteGaussianImageFilter);
78
80 using InputImageType = TInputImage;
81 using OutputImageType = TOutputImage;
82
85 using OutputPixelType = typename TOutputImage::PixelType;
86 using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
87 using InputPixelType = typename TInputImage::PixelType;
88 using InputInternalPixelType = typename TInputImage::InternalPixelType;
89
93
96 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
97
102
105#ifndef ITK_FUTURE_LEGACY_REMOVE
106 using InputBoundaryConditionPointerType [[deprecated("Please just use `BoundaryConditionType *` instead!")]] =
108#endif
112
116 using ScalarRealType = double;
117
121
129 itkSetMacro(Variance, ArrayType);
130 itkGetConstMacro(Variance, const ArrayType);
132
137 itkSetMacro(MaximumError, ArrayType);
138 itkGetConstMacro(MaximumError, const ArrayType);
140
144 itkGetConstMacro(MaximumKernelWidth, unsigned int);
145 itkSetMacro(MaximumKernelWidth, unsigned int);
147
154 itkGetConstMacro(FilterDimensionality, unsigned int);
155 itkSetMacro(FilterDimensionality, unsigned int);
157
160 itkSetMacro(InputBoundaryCondition, BoundaryConditionType *);
161 itkGetConstMacro(InputBoundaryCondition, BoundaryConditionType *);
162 itkSetMacro(RealBoundaryCondition, RealBoundaryConditionPointerType);
163 itkGetConstMacro(RealBoundaryCondition, RealBoundaryConditionPointerType);
165
168 void
170 {
171 m_Variance.Fill(v);
172 this->Modified();
173 }
174
175 void
177 {
178 m_MaximumError.Fill(v);
179 this->Modified();
180 }
181
182 void
183 SetVariance(const double * v)
184 {
185 ArrayType dv;
186
187 for (unsigned int i = 0; i < ImageDimension; ++i)
188 {
189 dv[i] = v[i];
190 }
191 this->SetVariance(dv);
192 }
193
194 void
195 SetVariance(const float * v)
196 {
197 ArrayType dv;
198
199 for (unsigned int i = 0; i < ImageDimension; ++i)
200 {
201 dv[i] = v[i];
202 }
203 this->SetVariance(dv);
204 }
205
208 void
209 SetSigma(const ArrayType & sigma)
210 {
211 ArrayType variance;
212 for (unsigned int i = 0; i < ImageDimension; ++i)
213 {
214 variance[i] = sigma[i] * sigma[i];
215 }
216 this->SetVariance(variance);
217 }
218
222 void
223 SetSigmaArray(const ArrayType & sigmas)
224 {
225 this->SetSigma(sigmas);
226 }
227 void
228 SetSigma(double sigma)
229 {
230 this->SetVariance(sigma * sigma);
231 }
232
233
235 ArrayType
237 {
238 ArrayType sigmas;
239 for (unsigned int i = 0; i < ImageDimension; ++i)
240 {
241 sigmas[i] = std::sqrt(m_Variance[i]);
242 }
243 return sigmas;
244 }
245
248 double
249 GetSigma() const
250 {
251 return std::sqrt(m_Variance[0]);
252 }
253
254 void
255 SetMaximumError(const double * v)
256 {
257 ArrayType dv;
258
259 for (unsigned int i = 0; i < ImageDimension; ++i)
260 {
261 dv[i] = v[i];
262 }
263 this->SetMaximumError(dv);
264 }
265
266 void
267 SetMaximumError(const float * v)
268 {
269 ArrayType dv;
270
271 for (unsigned int i = 0; i < ImageDimension; ++i)
272 {
273 dv[i] = v[i];
274 }
275 this->SetMaximumError(dv);
276 }
277
279 unsigned int
280 GetKernelRadius(const unsigned int dimension) const;
281
285
290
297 itkSetMacro(UseImageSpacing, bool);
298 itkGetConstMacro(UseImageSpacing, bool);
299 itkBooleanMacro(UseImageSpacing);
301
302#if !defined(ITK_FUTURE_LEGACY_REMOVE)
307 void
308 SetUseImageSpacingOn()
309 {
310 this->SetUseImageSpacing(true);
311 }
312
316 void
317 SetUseImageSpacingOff()
318 {
319 this->SetUseImageSpacing(false);
320 }
321#endif
322
332 itkLegacyMacro(unsigned int GetInternalNumberOfStreamDivisions() const;)
333 itkLegacyMacro(void SetInternalNumberOfStreamDivisions(unsigned int);)
334
336
337protected:
348
349 ~DiscreteGaussianImageFilter() override = default;
350 void
351 PrintSelf(std::ostream & os, Indent indent) const override;
352
359 void
361
367 void
368 GenerateData() override;
369
371 void
372 GenerateKernel(const unsigned int dimension, KernelType & oper) const;
373
377
378private:
382
387
390 unsigned int m_MaximumKernelWidth{};
391
394
397
401
404
407
410};
411} // end namespace itk
412
413#ifndef ITK_MANUAL_INSTANTIATION
414# include "itkDiscreteGaussianImageFilter.hxx"
415#endif
416
417#endif
void GenerateKernel(const unsigned int dimension, KernelType &oper) const
ZeroFluxNeumannBoundaryCondition< TInputImage > InputDefaultBoundaryConditionType
void GenerateInputRequestedRegion() override
typename NumericTraits< InputPixelType >::ValueType InputPixelValueType
typename TOutputImage::InternalPixelType OutputInternalPixelType
unsigned int GetKernelRadius(const unsigned int dimension) const
ImageToImageFilter< TInputImage, TOutputImage > Superclass
GaussianOperator< RealOutputPixelValueType, ImageDimension > KernelType
void SetVariance(const typename ArrayType::ValueType v)
Image< OutputPixelType, ImageDimension > RealOutputImageType
ImageBoundaryCondition< RealOutputImageType > * RealBoundaryConditionPointerType
ZeroFluxNeumannBoundaryCondition< RealOutputImageType > RealDefaultBoundaryConditionType
typename TInputImage::InternalPixelType InputInternalPixelType
unsigned int GetInternalNumberOfStreamDivisions() const
Set/Get number of pieces to divide the input for the internal composite pipeline. The upstream pipeli...
FixedArray< double, Self::ImageDimension > ArrayType
typename NumericTraits< RealOutputPixelType >::ValueType RealOutputPixelValueType
typename TInputImage::PixelType InputPixelType
void SetMaximumError(const typename ArrayType::ValueType v)
void SetInternalNumberOfStreamDivisions(unsigned int)
ArrayType GetKernelRadius() const
ArrayType GetKernelVarianceArray() const
~DiscreteGaussianImageFilter() override=default
typename NumericTraits< OutputPixelType >::ValueType OutputPixelValueType
typename TOutputImage::PixelType OutputPixelType
void PrintSelf(std::ostream &os, Indent indent) const override
ImageBoundaryCondition< TInputImage > BoundaryConditionType
typename NumericTraits< OutputPixelType >::RealType RealOutputPixelType
Simulate a standard C array with copy semantics.
A NeighborhoodOperator whose coefficients are a one dimensional, discrete Gaussian kernel.
A virtual base object that defines an interface to a class of boundary condition objects for use by n...
Templated n-dimensional image class.
Definition itkImage.h:89
Control indentation during Print() invocation.
Definition itkIndent.h:50
itk::Size< VDimension > RadiusType
virtual void Modified() const
Implements transparent reference counting.
A function object that determines a neighborhood of values at an image boundary according to a Neuman...
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....