1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285
|
/*=========================================================================
*
* Copyright Insight Software Consortium
*
* Licensed under the Apache License, Version 2.0 (the "License");
* you may not use this file except in compliance with the License.
* You may obtain a copy of the License at
*
* http://www.apache.org/licenses/LICENSE-2.0.txt
*
* Unless required by applicable law or agreed to in writing, software
* distributed under the License is distributed on an "AS IS" BASIS,
* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
* See the License for the specific language governing permissions and
* limitations under the License.
*
*=========================================================================*/
#ifndef itkDiscreteGaussianImageFilter_h
#define itkDiscreteGaussianImageFilter_h
#include "itkImageToImageFilter.h"
#include "itkImage.h"
namespace itk
{
/**
* \class DiscreteGaussianImageFilter
* \brief Blurs an image by separable convolution with discrete gaussian kernels.
* This filter performs Gaussian blurring by separable convolution of an image
* and a discrete Gaussian operator (kernel).
*
* The Gaussian operator used here was described by Tony Lindeberg (Discrete
* Scale-Space Theory and the Scale-Space Primal Sketch. Dissertation. Royal
* Institute of Technology, Stockholm, Sweden. May 1991.) The Gaussian kernel
* used here was designed so that smoothing and derivative operations commute
* after discretization.
*
* The variance or standard deviation (sigma) will be evaluated as pixel units
* if SetUseImageSpacing is off (false) or as physical units if
* SetUseImageSpacing is on (true, default). The variance can be set
* independently in each dimension.
*
* When the Gaussian kernel is small, this filter tends to run faster than
* itk::RecursiveGaussianImageFilter.
*
* \sa GaussianOperator
* \sa Image
* \sa Neighborhood
* \sa NeighborhoodOperator
* \sa RecursiveGaussianImageFilter
*
* \ingroup ImageEnhancement
* \ingroup ImageFeatureExtraction
* \ingroup ITKSmoothing
*
* \wiki
* \wikiexample{Smoothing/DiscreteGaussianImageFilter,Smooth an image with a discrete Gaussian filter}
* \endwiki
*/
template< typename TInputImage, typename TOutputImage >
class ITK_TEMPLATE_EXPORT DiscreteGaussianImageFilter:
public ImageToImageFilter< TInputImage, TOutputImage >
{
public:
/** Standard class typedefs. */
typedef DiscreteGaussianImageFilter Self;
typedef ImageToImageFilter< TInputImage, TOutputImage > Superclass;
typedef SmartPointer< Self > Pointer;
typedef SmartPointer< const Self > ConstPointer;
/** Method for creation through the object factory. */
itkNewMacro(Self);
/** Run-time type information (and related methods). */
itkTypeMacro(DiscreteGaussianImageFilter, ImageToImageFilter);
/** Image type information. */
typedef TInputImage InputImageType;
typedef TOutputImage OutputImageType;
/** Extract some information from the image types. Dimensionality
* of the two images is assumed to be the same. */
typedef typename TOutputImage::PixelType OutputPixelType;
typedef typename TOutputImage::InternalPixelType OutputInternalPixelType;
typedef typename TInputImage::PixelType InputPixelType;
typedef typename TInputImage::InternalPixelType InputInternalPixelType;
/** Pixel value type for Vector pixel types **/
typedef typename NumericTraits<InputPixelType>::ValueType InputPixelValueType;
typedef typename NumericTraits<OutputPixelType>::ValueType OutputPixelValueType;
/** Extract some information from the image types. Dimensionality
* of the two images is assumed to be the same. */
itkStaticConstMacro(ImageDimension, unsigned int,
TOutputImage::ImageDimension);
/** Typedef of double containers */
typedef FixedArray< double, itkGetStaticConstMacro(ImageDimension) > ArrayType;
/** The variance for the discrete Gaussian kernel. Sets the variance
* independently for each dimension, but
* see also SetVariance(const double v). The default is 0.0 in each
* dimension. If UseImageSpacing is true, the units are the physical units
* of your image. If UseImageSpacing is false then the units are
* pixels. */
itkSetMacro(Variance, ArrayType);
itkGetConstMacro(Variance, const ArrayType);
/** The algorithm will size the discrete kernel so that the error
* resulting from truncation of the kernel is no greater than
* MaximumError. The default is 0.01 in each dimension. */
itkSetMacro(MaximumError, ArrayType);
itkGetConstMacro(MaximumError, const ArrayType);
/** Set the kernel to be no wider than MaximumKernelWidth pixels,
* even if MaximumError demands it. The default is 32 pixels. */
itkGetConstMacro(MaximumKernelWidth, int);
itkSetMacro(MaximumKernelWidth, int);
/** Set the number of dimensions to smooth. Defaults to the image
* dimension. Can be set to less than ImageDimension, smoothing all
* the dimensions less than FilterDimensionality. For instance, to
* smooth the slices of a volume without smoothing in Z, set the
* FilterDimensionality to 2. */
itkGetConstMacro(FilterDimensionality, unsigned int);
itkSetMacro(FilterDimensionality, unsigned int);
/** Convenience Set methods for setting all dimensional parameters
* to the same values. */
void SetVariance(const typename ArrayType::ValueType v)
{
m_Variance.Fill(v);
this->Modified();
}
void SetMaximumError(const typename ArrayType::ValueType v)
{
m_MaximumError.Fill(v);
this->Modified();
}
void SetVariance(const double *v)
{
ArrayType dv;
for ( unsigned int i = 0; i < ImageDimension; i++ )
{
dv[i] = v[i];
}
this->SetVariance(dv);
}
void SetVariance(const float *v)
{
ArrayType dv;
for ( unsigned int i = 0; i < ImageDimension; i++ )
{
dv[i] = v[i];
}
this->SetVariance(dv);
}
void SetMaximumError(const double *v)
{
ArrayType dv;
for ( unsigned int i = 0; i < ImageDimension; i++ )
{
dv[i] = v[i];
}
this->SetMaximumError(dv);
}
void SetMaximumError(const float *v)
{
ArrayType dv;
for ( unsigned int i = 0; i < ImageDimension; i++ )
{
dv[i] = v[i];
}
this->SetMaximumError(dv);
}
/** Use the image spacing information in calculations. Use this option if you
* want to specify Gaussian variance in real world units. Default is
* ImageSpacingOn. */
void SetUseImageSpacingOn()
{ this->SetUseImageSpacing(true); }
/** Ignore the image spacing. Use this option if you want to specify Gaussian
variance in pixels. Default is ImageSpacingOn. */
void SetUseImageSpacingOff()
{ this->SetUseImageSpacing(false); }
/** Set/Get whether or not the filter will use the spacing of the input
image in its calculations */
itkSetMacro(UseImageSpacing, bool);
itkGetConstMacro(UseImageSpacing, bool);
/** \brief Set/Get number of pieces to divide the input for the
* internal composite pipeline. The upstream pipeline will not be
* effected.
*
* The default value is $ImageDimension^2$.
*
* This parameter was introduced to reduce the memory used by images
* internally, at the cost of performance.
*/
itkSetMacro(InternalNumberOfStreamDivisions, unsigned int);
itkGetConstReferenceMacro(InternalNumberOfStreamDivisions, unsigned int);
/** DiscreteGaussianImageFilter needs a larger input requested region
* than the output requested region (larger by the size of the
* Gaussian kernel). As such, DiscreteGaussianImageFilter needs to
* provide an implementation for GenerateInputRequestedRegion() in
* order to inform the pipeline execution model.
* \sa ImageToImageFilter::GenerateInputRequestedRegion() */
virtual void GenerateInputRequestedRegion() ITK_OVERRIDE;
#ifdef ITK_USE_CONCEPT_CHECKING
// Begin concept checking
itkConceptMacro( OutputHasNumericTraitsCheck,
( Concept::HasNumericTraits< OutputPixelValueType > ) );
// End concept checking
#endif
protected:
DiscreteGaussianImageFilter()
{
m_Variance.Fill(0.0);
m_MaximumError.Fill(0.01);
m_MaximumKernelWidth = 32;
m_UseImageSpacing = true;
m_FilterDimensionality = ImageDimension;
m_InternalNumberOfStreamDivisions = ImageDimension * ImageDimension;
}
virtual ~DiscreteGaussianImageFilter() ITK_OVERRIDE {}
void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
/** Standard pipeline method. While this class does not implement a
* ThreadedGenerateData(), its GenerateData() delegates all
* calculations to an NeighborhoodOperatorImageFilter. Since the
* NeighborhoodOperatorImageFilter is multithreaded, this filter is
* multithreaded by default. */
void GenerateData() ITK_OVERRIDE;
private:
ITK_DISALLOW_COPY_AND_ASSIGN(DiscreteGaussianImageFilter);
/** The variance of the gaussian blurring kernel in each dimensional
direction. */
ArrayType m_Variance;
/** The maximum error of the gaussian blurring kernel in each dimensional
* direction. For definition of maximum error, see GaussianOperator.
* \sa GaussianOperator */
ArrayType m_MaximumError;
/** Maximum allowed kernel width for any dimension of the discrete Gaussian
approximation */
int m_MaximumKernelWidth;
/** Number of dimensions to process. Default is all dimensions */
unsigned int m_FilterDimensionality;
/** Flag to indicate whether to use image spacing */
bool m_UseImageSpacing;
/** Number of pieces to divide the input on the internal composite
pipeline. The upstream pipeline will not be effected. */
unsigned int m_InternalNumberOfStreamDivisions;
};
} // end namespace itk
#ifndef ITK_MANUAL_INSTANTIATION
#include "itkDiscreteGaussianImageFilter.hxx"
#endif
#endif
|