ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkConfidenceConnectedImageFilter.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 itkConfidenceConnectedImageFilter_h
19 #define itkConfidenceConnectedImageFilter_h
20 
21 #include "itkImage.h"
22 #include "itkImageToImageFilter.h"
23 
24 namespace itk
25 {
62 template< typename TInputImage, typename TOutputImage >
63 class ITK_TEMPLATE_EXPORT ConfidenceConnectedImageFilter:
64  public ImageToImageFilter< TInputImage, TOutputImage >
65 {
66 public:
72 
74  itkNewMacro(Self);
75 
77  itkTypeMacro(ConfidenceConnectedImageFilter,
79 
80  typedef TInputImage InputImageType;
81  typedef typename InputImageType::Pointer InputImagePointer;
82  typedef typename InputImageType::RegionType InputImageRegionType;
83  typedef typename InputImageType::PixelType InputImagePixelType;
86 
87  typedef TOutputImage OutputImageType;
88  typedef typename OutputImageType::Pointer OutputImagePointer;
89  typedef typename OutputImageType::RegionType OutputImageRegionType;
90  typedef typename OutputImageType::PixelType OutputImagePixelType;
91 
92  typedef std::vector< IndexType > SeedsContainerType;
93 
94  typedef typename NumericTraits<
96 
97  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
98 
100  void SetSeed(const IndexType & seed);
101 
103  void ClearSeeds();
104 
106  void AddSeed(const IndexType & seed);
107 
110  itkSetMacro(Multiplier, double);
111  itkGetConstMacro(Multiplier, double);
113 
115  itkSetMacro(NumberOfIterations, unsigned int);
116  itkGetConstMacro(NumberOfIterations, unsigned int);
118 
120  itkSetMacro(ReplaceValue, OutputImagePixelType);
121  itkGetConstMacro(ReplaceValue, OutputImagePixelType);
123 
126  itkSetMacro(InitialNeighborhoodRadius, unsigned int);
127  itkGetConstReferenceMacro(InitialNeighborhoodRadius, unsigned int);
129 
133  itkGetConstReferenceMacro(Mean, InputRealType);
134 
138  itkGetConstReferenceMacro(Variance, InputRealType);
139 
141  virtual const SeedsContainerType &GetSeeds() const;
142 
143 #ifdef ITK_USE_CONCEPT_CHECKING
144  // Begin concept checking
145  itkConceptMacro( InputHasNumericTraitsCheck,
147  itkConceptMacro( OutputHasNumericTraitsCheck,
149  // End concept checking
150 #endif
151 
152 protected:
155 
156  // Override since the filter needs all the data for the algorithm
157  void GenerateInputRequestedRegion() ITK_OVERRIDE;
158 
159  // Override since the filter produces the entire dataset
160  void EnlargeOutputRequestedRegion(DataObject *output) ITK_OVERRIDE;
161 
162  void GenerateData() ITK_OVERRIDE;
163 
164 private:
165  ITK_DISALLOW_COPY_AND_ASSIGN(ConfidenceConnectedImageFilter);
166 
167  SeedsContainerType m_Seeds;
168  double m_Multiplier;
169  unsigned int m_NumberOfIterations;
170  OutputImagePixelType m_ReplaceValue;
171  unsigned int m_InitialNeighborhoodRadius;
173  InputRealType m_Variance;
174 };
175 } // end namespace itk
176 
177 #ifndef ITK_MANUAL_INSTANTIATION
178 #include "itkConfidenceConnectedImageFilter.hxx"
179 #endif
180 
181 #endif
NumericTraits< InputImagePixelType >::RealType InputRealType
ImageToImageFilter< TInputImage, TOutputImage > Superclass
Base class for all process objects that output image data.
OutputImageType::PixelType OutputImagePixelType
Segment pixels with similar statistics using connectivity.
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
Define additional traits for native types such as int or float.
#define itkConceptMacro(name, concept)
Base class for all data objects in ITK.