ITK  6.0.0
Insight Toolkit
itkScalarImageToCooccurrenceMatrixFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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  * https://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  typename TMaskImageType = TImageType>
96 class ITK_TEMPLATE_EXPORT ScalarImageToCooccurrenceMatrixFilter : public ProcessObject
97 {
98 public:
99  ITK_DISALLOW_COPY_AND_MOVE(ScalarImageToCooccurrenceMatrixFilter);
100 
106 
108  itkOverrideGetNameOfClassMacro(ScalarImageToCooccurrenceMatrixFilter);
109 
111  itkNewMacro(Self);
112 
113  using ImageType = TImageType;
116  using PixelType = typename ImageType::PixelType;
118  using RadiusType = typename ImageType::SizeType;
119  using OffsetType = typename ImageType::OffsetType;
123  using MaskImageType = TMaskImageType;
126  using MaskPixelType = typename MaskImageType::PixelType;
127 
129 
134 
135  static constexpr unsigned int DefaultBinsPerAxis = 256;
136 
139  itkSetConstObjectMacro(Offsets, OffsetVector);
140  itkGetConstObjectMacro(Offsets, OffsetVector);
143  void
144  SetOffset(const OffsetType offset);
145 
147  itkSetMacro(NumberOfBinsPerAxis, unsigned int);
148  itkGetConstMacro(NumberOfBinsPerAxis, unsigned int);
153  void
154  SetPixelValueMinMax(PixelType min, PixelType max);
155 
156  itkGetConstMacro(Min, PixelType);
157  itkGetConstMacro(Max, PixelType);
158 
161  itkSetMacro(Normalize, bool);
162  itkGetConstMacro(Normalize, bool);
163  itkBooleanMacro(Normalize);
167  using Superclass::SetInput;
168  void
169  SetInput(const ImageType * image);
170 
171  const ImageType *
172  GetInput() const;
173 
175  void
176  SetMaskImage(const MaskImageType * image);
177 
178  const MaskImageType *
179  GetMaskImage() const;
180 
182  const HistogramType *
183  GetOutput() const;
184 
187  itkSetMacro(InsidePixelValue, MaskPixelType);
188  itkGetConstMacro(InsidePixelValue, MaskPixelType);
191 protected:
193  ~ScalarImageToCooccurrenceMatrixFilter() override = default;
194  void
195  PrintSelf(std::ostream & os, Indent indent) const override;
196 
197  virtual void
198  FillHistogram(RadiusType radius, RegionType region);
199 
200  virtual void
201  FillHistogramWithMask(RadiusType radius, RegionType region, const MaskImageType * maskImage);
202 
205 
207  using Superclass::MakeOutput;
209  MakeOutput(DataObjectPointerArraySizeType idx) override;
210 
212  void
213  GenerateData() override;
214 
215 private:
216  void
217  NormalizeHistogram();
218 
220  PixelType m_Min{};
221  PixelType m_Max{};
222 
223  unsigned int m_NumberOfBinsPerAxis{};
224  MeasurementVectorType m_LowerBound{};
225  MeasurementVectorType m_UpperBound{};
226  bool m_Normalize{};
227 
228  MaskPixelType m_InsidePixelValue{};
229 };
230 } // end of namespace Statistics
231 } // end of namespace itk
232 
233 #ifndef ITK_MANUAL_INSTANTIATION
234 # include "itkScalarImageToCooccurrenceMatrixFilter.hxx"
235 #endif
236 
237 #endif
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::OffsetType
typename ImageType::OffsetType OffsetType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:119
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::PixelType
typename ImageType::PixelType PixelType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:116
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::MaskImageType
TMaskImageType MaskImageType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:123
ConstPointer
SmartPointer< const Self > ConstPointer
Definition: itkAddImageFilter.h:94
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::MaskPointer
typename MaskImageType::Pointer MaskPointer
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:124
itk::detail::VectorContainer
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
Definition: itkVectorContainer.h:51
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::OffsetVectorConstPointer
typename OffsetVector::ConstPointer OffsetVectorConstPointer
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:122
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itkImage.h
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ProcessObject::DataObjectPointerArraySizeType
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Definition: itkProcessObject.h:194
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::RadiusType
typename ImageType::SizeType RadiusType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:118
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::MaskPixelType
typename MaskImageType::PixelType MaskPixelType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:126
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter
This class computes a co-occurrence matrix (histogram) from a given image and a mask image if provide...
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:96
itk::Statistics::Histogram
This class stores measurement vectors in the context of n-dimensional histogram.
Definition: itkHistogram.h:77
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:55
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::HistogramPointer
typename HistogramType::Pointer HistogramPointer
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:131
itkHistogram.h
itkProcessObject.h
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::MeasurementType
typename NumericTraits< PixelType >::RealType MeasurementType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:128
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::MeasurementVectorType
typename HistogramType::MeasurementVectorType MeasurementVectorType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:133
itk::NumericTraits
Define additional traits for native types such as int or float.
Definition: itkNumericTraits.h:60
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::HistogramConstPointer
typename HistogramType::ConstPointer HistogramConstPointer
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:132
itkVectorContainer.h
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:139
itk::ProcessObject
class ITK_FORWARD_EXPORT ProcessObject
Definition: itkDataObject.h:41
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::OffsetVectorPointer
typename OffsetVector::Pointer OffsetVectorPointer
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:121
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::ImagePointer
typename ImageType::Pointer ImagePointer
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:114
itk::Array
Array class with size defined at construction time.
Definition: itkArray.h:47
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::RegionType
typename ImageType::RegionType RegionType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:117
itkNumericTraits.h
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::MaskConstPointer
typename MaskImageType::ConstPointer MaskConstPointer
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:125
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::ImageConstPointer
typename ImageType::ConstPointer ImageConstPointer
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:115
itk::Statistics::ScalarImageToCooccurrenceMatrixFilter::ImageType
TImageType ImageType
Definition: itkScalarImageToCooccurrenceMatrixFilter.h:113
itk::DataObject::Pointer
SmartPointer< Self > Pointer
Definition: itkDataObject.h:301