ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkNoiseImageFilter.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 itkNoiseImageFilter_h
19 #define itkNoiseImageFilter_h
20 
21 #include "itkBoxImageFilter.h"
22 #include "itkImage.h"
23 #include "itkNumericTraits.h"
24 
25 namespace itk
26 {
50 template< typename TInputImage, typename TOutputImage >
52  public BoxImageFilter< TInputImage, TOutputImage >
53 {
54 public:
56  itkStaticConstMacro(InputImageDimension, unsigned int,
57  TInputImage::ImageDimension);
58  itkStaticConstMacro(OutputImageDimension, unsigned int,
59  TOutputImage::ImageDimension);
61 
63  typedef TInputImage InputImageType;
64  typedef TOutputImage OutputImageType;
65 
71 
73  itkNewMacro(Self);
74 
77 
79  typedef typename InputImageType::PixelType InputPixelType;
80  typedef typename OutputImageType::PixelType OutputPixelType;
82 
83  typedef typename InputImageType::RegionType InputImageRegionType;
84  typedef typename OutputImageType::RegionType OutputImageRegionType;
85 
86  typedef typename InputImageType::SizeType InputSizeType;
87 
88 #ifdef ITK_USE_CONCEPT_CHECKING
89  // Begin concept checking
90  itkConceptMacro( InputHasNumericTraitsCheck,
92  // End concept checking
93 #endif
94 
95 protected:
97  virtual ~NoiseImageFilter() {}
98 
109  void ThreadedGenerateData(const OutputImageRegionType & outputRegionForThread,
110  ThreadIdType threadId) ITK_OVERRIDE;
111 
112 private:
113  NoiseImageFilter(const Self &); //purposely not implemented
114  void operator=(const Self &); //purposely not implemented
115 };
116 } // end namespace itk
117 
118 #ifndef ITK_MANUAL_INSTANTIATION
119 #include "itkNoiseImageFilter.hxx"
120 #endif
121 
122 #endif
InputImageType::SizeType InputSizeType
void operator=(const Self &)
InputImageType::RegionType InputImageRegionType
Base class for all process objects that output image data.
OutputImageType::RegionType OutputImageRegionType
InputImageType::PixelType InputPixelType
SmartPointer< Self > Pointer
BoxImageFilter< InputImageType, OutputImageType > Superclass
A base class for all the filters working on a box neighborhood.
NumericTraits< InputPixelType >::RealType InputRealType
OutputImageType::PixelType OutputPixelType
unsigned int ThreadIdType
Definition: itkIntTypes.h:159
SmartPointer< const Self > ConstPointer
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, ThreadIdType threadId) override
static const unsigned int OutputImageDimension
Define additional traits for native types such as int or float.
Calculate the local noise in an image.
#define itkConceptMacro(name, concept)
static const unsigned int InputImageDimension