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
27namespace itk
28{
29namespace Statistics
30{
93template <typename TImageType,
94 typename THistogramFrequencyContainer = DenseFrequencyContainer2,
95 typename TMaskImageType = TImageType>
96class ITK_TEMPLATE_EXPORT ScalarImageToCooccurrenceMatrixFilter : public ProcessObject
97{
98public:
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;
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
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
177
178 const MaskImageType *
180
182 const HistogramType *
183 GetOutput() const;
184
187 itkSetMacro(InsidePixelValue, MaskPixelType);
188 itkGetConstMacro(InsidePixelValue, MaskPixelType);
191protected:
194 void
195 PrintSelf(std::ostream & os, Indent indent) const override;
196
197 virtual void
199
200 virtual void
201 FillHistogramWithMask(RadiusType radius, RegionType region, const MaskImageType * maskImage);
202
205
207 using Superclass::MakeOutput;
210
212 void
213 GenerateData() override;
214
215private:
216 void
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
SmartPointer< Self > Pointer
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
Define additional traits for native types such as int or float.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
This class stores measurement vectors in the context of n-dimensional histogram.
Definition: itkHistogram.h:78
This class computes a co-occurrence matrix (histogram) from a given image and a mask image if provide...
virtual void FillHistogram(RadiusType radius, RegionType region)
void SetPixelValueMinMax(PixelType min, PixelType max)
void PrintSelf(std::ostream &os, Indent indent) const override
virtual void FillHistogramWithMask(RadiusType radius, RegionType region, const MaskImageType *maskImage)
DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override
Make a DataObject of the correct type to used as the specified output.
Define a front-end to the STL "vector" container that conforms to the IndexedContainerInterface.
SmartPointer< const Self > ConstPointer
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
class ITK_FORWARD_EXPORT ProcessObject
Definition: itkDataObject.h:41