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
128 itkSetMacro(Variance, ArrayType);
129 itkGetConstMacro(Variance, const ArrayType);
131
135 itkSetMacro(MaximumError, ArrayType);
136 itkGetConstMacro(MaximumError, const ArrayType);
138
141 itkGetConstMacro(MaximumKernelWidth, unsigned int);
142 itkSetMacro(MaximumKernelWidth, unsigned int);
144
150 itkGetConstMacro(FilterDimensionality, unsigned int);
151 itkSetMacro(FilterDimensionality, unsigned int);
153
155 itkSetMacro(InputBoundaryCondition, BoundaryConditionType *);
156 itkGetConstMacro(InputBoundaryCondition, BoundaryConditionType *);
157 itkSetMacro(RealBoundaryCondition, RealBoundaryConditionPointerType);
158 itkGetConstMacro(RealBoundaryCondition, RealBoundaryConditionPointerType);
160
163 void
165 {
166 m_Variance.Fill(v);
167 this->Modified();
168 }
169
170
171 void
173 {
174 m_MaximumError.Fill(v);
175 this->Modified();
176 }
177
178 void
179 SetVariance(const double * v)
180 {
181 ArrayType dv;
182
183 for (unsigned int i = 0; i < ImageDimension; ++i)
184 {
185 dv[i] = v[i];
186 }
187 this->SetVariance(dv);
188 }
189
190 void
191 SetVariance(const float * v)
192 {
193 ArrayType dv;
194
195 for (unsigned int i = 0; i < ImageDimension; ++i)
196 {
197 dv[i] = v[i];
198 }
199 this->SetVariance(dv);
200 }
201
204 void
205 SetSigma(const ArrayType & sigma)
206 {
207 ArrayType variance;
208 for (unsigned int i = 0; i < ImageDimension; ++i)
209 {
210 variance[i] = sigma[i] * sigma[i];
211 }
212 this->SetVariance(variance);
213 }
214
217 void
218 SetSigmaArray(const ArrayType & sigmas)
219 {
220 this->SetSigma(sigmas);
221 }
222 void
223 SetSigma(double sigma)
224 {
225 this->SetVariance(sigma * sigma);
226 }
227
228
230 ArrayType
232 {
233 ArrayType sigmas;
234 for (unsigned int i = 0; i < ImageDimension; ++i)
235 {
236 sigmas[i] = std::sqrt(m_Variance[i]);
237 }
238 return sigmas;
239 }
240
241
244 double
245 GetSigma() const
246 {
247 return std::sqrt(m_Variance[0]);
248 }
249
250 void
251 SetMaximumError(const double * v)
252 {
253 ArrayType dv;
254
255 for (unsigned int i = 0; i < ImageDimension; ++i)
256 {
257 dv[i] = v[i];
258 }
259 this->SetMaximumError(dv);
260 }
261
262 void
263 SetMaximumError(const float * v)
264 {
265 ArrayType dv;
266
267 for (unsigned int i = 0; i < ImageDimension; ++i)
268 {
269 dv[i] = v[i];
270 }
271 this->SetMaximumError(dv);
272 }
273
275 unsigned int
276 GetKernelRadius(const unsigned int dimension) const;
277
281
286
292 itkSetMacro(UseImageSpacing, bool);
293 itkGetConstMacro(UseImageSpacing, bool);
294 itkBooleanMacro(UseImageSpacing);
296
297#if !defined(ITK_FUTURE_LEGACY_REMOVE)
302 void
303 SetUseImageSpacingOn()
304 {
305 this->SetUseImageSpacing(true);
306 }
307
311 void
312 SetUseImageSpacingOff()
313 {
314 this->SetUseImageSpacing(false);
315 }
316#endif
317
327 itkLegacyMacro(unsigned int GetInternalNumberOfStreamDivisions() const;)
328 itkLegacyMacro(void SetInternalNumberOfStreamDivisions(unsigned int);)
329
331
332protected:
343
344 ~DiscreteGaussianImageFilter() override = default;
345 void
346 PrintSelf(std::ostream & os, Indent indent) const override;
347
354 void
356
362 void
363 GenerateData() override;
364
366 void
367 GenerateKernel(const unsigned int dimension, KernelType & oper) const;
368
372
373private:
377
382
385 unsigned int m_MaximumKernelWidth{};
386
389
392
396
399
402
405};
406} // end namespace itk
407
408#ifndef ITK_MANUAL_INSTANTIATION
409# include "itkDiscreteGaussianImageFilter.hxx"
410#endif
411
412#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....