ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkScalarImageToCooccurrenceMatrixFilter.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 itkScalarImageToCooccurrenceMatrixFilter_h
19 #define itkScalarImageToCooccurrenceMatrixFilter_h
20 
21 #include "itkImage.h"
22 #include "itkHistogram.h"
23 #include "itkVectorContainer.h"
24 #include "itkNumericTraits.h"
25 #include "itkProcessObject.h"
26 
27 namespace itk
28 {
29 namespace Statistics
30 {
93 template< typename TImageType,
94  typename THistogramFrequencyContainer = DenseFrequencyContainer2 >
95 class ITK_TEMPLATE_EXPORT ScalarImageToCooccurrenceMatrixFilter:public ProcessObject
96 {
97 public:
103 
106 
108  itkNewMacro(Self);
109 
110  typedef TImageType ImageType;
111  typedef typename ImageType::Pointer ImagePointer;
112  typedef typename ImageType::ConstPointer ImageConstPointer;
113  typedef typename ImageType::PixelType PixelType;
114  typedef typename ImageType::RegionType RegionType;
115  typedef typename ImageType::SizeType RadiusType;
116  typedef typename ImageType::OffsetType OffsetType;
120 
122 
127 
128  itkStaticConstMacro(DefaultBinsPerAxis, unsigned int, 256);
129 
132  itkSetConstObjectMacro(Offsets, OffsetVector);
133  itkGetConstObjectMacro(Offsets, OffsetVector);
135 
136  void SetOffset(const OffsetType offset);
137 
139  itkSetMacro(NumberOfBinsPerAxis, unsigned int);
140  itkGetConstMacro(NumberOfBinsPerAxis, unsigned int);
142 
145  void SetPixelValueMinMax(PixelType min, PixelType max);
146 
147  itkGetConstMacro(Min, PixelType);
148  itkGetConstMacro(Max, PixelType);
149 
152  itkSetMacro(Normalize, bool);
153  itkGetConstMacro(Normalize, bool);
154  itkBooleanMacro(Normalize);
156 
158  using Superclass::SetInput;
159  void SetInput(const ImageType *image);
160 
161  const ImageType * GetInput() const;
162 
164  void SetMaskImage(const ImageType *image);
165 
166  const ImageType * GetMaskImage() const;
167 
169  const HistogramType * GetOutput() const;
170 
173  itkSetMacro(InsidePixelValue, PixelType);
174  itkGetConstMacro(InsidePixelValue, PixelType);
176 
177 protected:
179  virtual ~ScalarImageToCooccurrenceMatrixFilter() ITK_OVERRIDE {}
180  virtual void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
181 
182  virtual void FillHistogram(RadiusType radius, RegionType region);
183 
184  virtual void FillHistogramWithMask(RadiusType radius, RegionType region, const ImageType *maskImage);
185 
188 
190  using Superclass::MakeOutput;
191  virtual DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) ITK_OVERRIDE;
192 
194  virtual void GenerateData() ITK_OVERRIDE;
195 
196 private:
197  ITK_DISALLOW_COPY_AND_ASSIGN(ScalarImageToCooccurrenceMatrixFilter);
198 
199  void NormalizeHistogram();
200 
202  PixelType m_Min;
203  PixelType m_Max;
204 
205  unsigned int m_NumberOfBinsPerAxis;
206  MeasurementVectorType m_LowerBound;
207  MeasurementVectorType m_UpperBound;
208  bool m_Normalize;
209 
210  PixelType m_InsidePixelValue;
211 };
212 } // end of namespace Statistics
213 } // end of namespace itk
214 
215 #ifndef ITK_MANUAL_INSTANTIATION
216 #include "itkScalarImageToCooccurrenceMatrixFilter.hxx"
217 #endif
218 
219 #endif
This class computes a co-occurence matrix (histogram) from a given image and a mask image if provided...
Light weight base class for most itk classes.
This class stores measurement vectors in the context of n-dimensional histogram.
Definition: itkHistogram.h:77
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Histogram< MeasurementType, THistogramFrequencyContainer > HistogramType
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
Superclass::MeasurementVectorType MeasurementVectorType
Definition: itkHistogram.h:101
Define a front-end to the STL &quot;vector&quot; container that conforms to the IndexedContainerInterface.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Define additional traits for native types such as int or float.