ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkDiscreteGaussianDerivativeImageFilter.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 itkDiscreteGaussianDerivativeImageFilter_h
19 #define itkDiscreteGaussianDerivativeImageFilter_h
20 
21 #include "itkImageToImageFilter.h"
22 #include "itkImage.h"
23 
24 namespace itk
25 {
58 template< typename TInputImage, typename TOutputImage >
59 class ITK_TEMPLATE_EXPORT DiscreteGaussianDerivativeImageFilter :
60  public ImageToImageFilter< TInputImage, TOutputImage >
61 {
62 public:
63  ITK_DISALLOW_COPY_AND_ASSIGN(DiscreteGaussianDerivativeImageFilter);
64 
70 
72  itkNewMacro(Self);
73 
76 
78  using InputImageType = TInputImage;
79  using OutputImageType = TOutputImage;
80 
83  using OutputPixelType = typename TOutputImage::PixelType;
84  using OutputInternalPixelType = typename TOutputImage::InternalPixelType;
85  using InputPixelType = typename TInputImage::PixelType;
86  using InputInternalPixelType = typename TInputImage::InternalPixelType;
87 
90  static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
91 
94 
97 
101  itkSetMacro(Order, OrderArrayType);
102  itkGetConstMacro(Order, const OrderArrayType);
104 
111  itkSetMacro(Variance, ArrayType);
112  itkGetConstMacro(Variance, const ArrayType);
114 
118  itkSetMacro(MaximumError, ArrayType);
119  itkGetConstMacro(MaximumError, const ArrayType);
121 
124  itkGetConstMacro(MaximumKernelWidth, int);
125  itkSetMacro(MaximumKernelWidth, int);
127 
137  itkSetMacro(InternalNumberOfStreamDivisions, unsigned int);
138  itkGetConstMacro(InternalNumberOfStreamDivisions, unsigned int);
139 
144  void SetOrder(const typename OrderArrayType::ValueType v)
145  {
146  OrderArrayType a;
147 
148  a.Fill(v);
149  this->SetOrder(a);
150  }
151 
152  void SetVariance(const typename ArrayType::ValueType v)
153  {
154  ArrayType a;
155 
156  a.Fill(v);
157  this->SetVariance(a);
158  }
159 
160  void SetMaximumError(const typename ArrayType::ValueType v)
161  {
162  ArrayType a;
163 
164  a.Fill(v);
165  this->SetMaximumError(a);
166  }
167 
172  itkSetMacro(UseImageSpacing, bool);
173  itkGetConstMacro(UseImageSpacing, bool);
174  itkBooleanMacro(UseImageSpacing);
176 
180  itkSetMacro(NormalizeAcrossScale, bool);
181  itkGetConstMacro(NormalizeAcrossScale, bool);
182  itkBooleanMacro(NormalizeAcrossScale);
184 
185 #ifdef ITK_USE_CONCEPT_CHECKING
186  // Begin concept checking
187  itkConceptMacro( OutputHasNumericTraitsCheck,
189  // End concept checking
190 #endif
191 
192 protected:
193 
195  {
196  m_Order.Fill(1);
197  m_Variance.Fill(0.0);
198  m_MaximumError.Fill(0.01);
199  m_MaximumKernelWidth = 32;
200  m_UseImageSpacing = true;
201  m_NormalizeAcrossScale = false;
202  m_InternalNumberOfStreamDivisions = ImageDimension * ImageDimension;
203  }
204 
205  ~DiscreteGaussianDerivativeImageFilter() override = default;
206  void
207  PrintSelf(std::ostream & os, Indent indent) const override;
208 
215  void
216  GenerateInputRequestedRegion() override;
217 
223  void
224  GenerateData() override;
225 
226 private:
227 
230 
234 
239 
243 
246 
249 
253 };
254 } // end namespace itk
255 
256 #ifndef ITK_MANUAL_INSTANTIATION
257 #include "itkDiscreteGaussianDerivativeImageFilter.hxx"
258 #endif
259 
260 #endif
void SetOrder(const typename OrderArrayType::ValueType v)
typename TOutputImage::InternalPixelType OutputInternalPixelType
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.
void Fill(const ValueType &)
TOutputImage OutputImageType
void SetMaximumError(const typename ArrayType::ValueType v)
typename TInputImage::InternalPixelType InputInternalPixelType
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
Calculates image derivatives using discrete derivative gaussian kernels. This filter calculates Gauss...
#define itkConceptMacro(name, concept)