ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkDiscreteGaussianImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
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  * http://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 "itkImageToImageFilter.h"
22 #include "itkImage.h"
23 
24 namespace itk
25 {
61 template< typename TInputImage, typename TOutputImage >
62 class ITK_TEMPLATE_EXPORT DiscreteGaussianImageFilter:
63  public ImageToImageFilter< TInputImage, TOutputImage >
64 {
65 public:
66  ITK_DISALLOW_COPY_AND_ASSIGN(DiscreteGaussianImageFilter);
67 
73 
75  itkNewMacro(Self);
76 
79 
81  using InputImageType = TInputImage;
82  using OutputImageType = TOutputImage;
83 
86  using OutputPixelType = typename TOutputImage::PixelType;
87  using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
88  using InputPixelType = typename TInputImage::PixelType;
89  using InputInternalPixelType = typename TInputImage::InternalPixelType;
90 
94 
97  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
98 
101 
108  itkSetMacro(Variance, ArrayType);
109  itkGetConstMacro(Variance, const ArrayType);
111 
115  itkSetMacro(MaximumError, ArrayType);
116  itkGetConstMacro(MaximumError, const ArrayType);
118 
121  itkGetConstMacro(MaximumKernelWidth, int);
122  itkSetMacro(MaximumKernelWidth, int);
124 
130  itkGetConstMacro(FilterDimensionality, unsigned int);
131  itkSetMacro(FilterDimensionality, unsigned int);
133 
136  void SetVariance(const typename ArrayType::ValueType v)
137  {
138  m_Variance.Fill(v);
139  this->Modified();
140  }
142 
143  void SetMaximumError(const typename ArrayType::ValueType v)
144  {
145  m_MaximumError.Fill(v);
146  this->Modified();
147  }
148 
149  void SetVariance(const double *v)
150  {
151  ArrayType dv;
152 
153  for ( unsigned int i = 0; i < ImageDimension; i++ )
154  {
155  dv[i] = v[i];
156  }
157  this->SetVariance(dv);
158  }
159 
160  void SetVariance(const float *v)
161  {
162  ArrayType dv;
163 
164  for ( unsigned int i = 0; i < ImageDimension; i++ )
165  {
166  dv[i] = v[i];
167  }
168  this->SetVariance(dv);
169  }
170 
171  void SetMaximumError(const double *v)
172  {
173  ArrayType dv;
174 
175  for ( unsigned int i = 0; i < ImageDimension; i++ )
176  {
177  dv[i] = v[i];
178  }
179  this->SetMaximumError(dv);
180  }
181 
182  void SetMaximumError(const float *v)
183  {
184  ArrayType dv;
185 
186  for ( unsigned int i = 0; i < ImageDimension; i++ )
187  {
188  dv[i] = v[i];
189  }
190  this->SetMaximumError(dv);
191  }
192 
197  { this->SetUseImageSpacing(true); }
198 
202  { this->SetUseImageSpacing(false); }
203 
206  itkSetMacro(UseImageSpacing, bool);
207  itkGetConstMacro(UseImageSpacing, bool);
209 
219  itkLegacyMacro(unsigned int GetInternalNumberOfStreamDivisions() const);
220  itkLegacyMacro(void SetInternalNumberOfStreamDivisions(unsigned int));
221 
228  void GenerateInputRequestedRegion() override;
229 
230 #ifdef ITK_USE_CONCEPT_CHECKING
231  // Begin concept checking
232 
233  itkConceptMacro( OutputHasNumericTraitsCheck,
235 
236  // End concept checking
237 #endif
238 
239 protected:
241  {
242  m_Variance.Fill(0.0);
243  m_MaximumError.Fill(0.01);
244  m_MaximumKernelWidth = 32;
245  m_UseImageSpacing = true;
246  m_FilterDimensionality = ImageDimension;
247  }
248 
249  ~DiscreteGaussianImageFilter() override = default;
250  void PrintSelf(std::ostream & os, Indent indent) const override;
251 
257  void GenerateData() override;
258 
259 private:
263 
268 
272 
275 
278 
279 };
280 } // end namespace itk
281 
282 #ifndef ITK_MANUAL_INSTANTIATION
283 #include "itkDiscreteGaussianImageFilter.hxx"
284 #endif
285 
286 #endif
typename TInputImage::InternalPixelType InputInternalPixelType
Define numeric traits for std::vector.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Base class for all process objects that output image data.
typename TOutputImage::PixelType OutputPixelType
void SetVariance(const typename ArrayType::ValueType v)
typename TInputImage::PixelType InputPixelType
typename NumericTraits< InputPixelType >::ValueType InputPixelValueType
typename NumericTraits< OutputPixelType >::ValueType OutputPixelValueType
TOutputImage OutputImageType
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
typename TOutputImage::InternalPixelType OutputInternalPixelType
#define itkConceptMacro(name, concept)
Blurs an image by separable convolution with discrete gaussian kernels. This filter performs Gaussian...
void SetMaximumError(const typename ArrayType::ValueType v)