ITK  6.0.0
Insight Toolkit
itkScalarImageToCooccurrenceListSampleFilter.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 itkScalarImageToCooccurrenceListSampleFilter_h
19#define itkScalarImageToCooccurrenceListSampleFilter_h
20
21#include <typeinfo>
22
24#include "itkSmartPointer.h"
29#include "itkFixedArray.h"
30
31#include <vector>
32#include <algorithm>
33#include <iostream>
34
35namespace itk
36{
37namespace Statistics
38{
51template <typename TImage>
53{
54public:
55 ITK_DISALLOW_COPY_AND_MOVE(ScalarImageToCooccurrenceListSampleFilter);
56
57 using ImageType = TImage;
58
60
63
69
72
75 using OffsetTable = std::vector<OffsetType>;
76
77 void
78 UseNeighbor(const OffsetType & offset);
79
81 using Superclass::SetInput;
82 void
83 SetInput(const ImageType * image);
84
85 const ImageType *
86 GetInput() const;
87
89 const SampleType *
90 GetOutput() const;
91
93 itkOverrideGetNameOfClassMacro(ScalarImageToCooccurrenceListSampleFilter);
94
96 itkNewMacro(Self);
97
99 static constexpr unsigned int MeasurementVectorSize = 2;
100
102 static constexpr unsigned int ImageDimension = TImage::ImageDimension;
103
104protected:
107 void
108 PrintSelf(std::ostream & os, Indent indent) const override;
109
112 using Superclass::MakeOutput;
115
117 void
118 GenerateData() override;
119
120private:
121 OffsetTable m_OffsetTable{};
122}; // end of class ScalarImageToListSampleFilter
123} // end of namespace Statistics
124} // end of namespace itk
125
126#ifndef ITK_MANUAL_INSTANTIATION
127# include "itkScalarImageToCooccurrenceListSampleFilter.hxx"
128#endif
129
130#endif
SmartPointer< Self > Pointer
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:54
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Light weight base class for most itk classes.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
A neighborhood iterator which can take on an arbitrary shape.
This class is the native implementation of the a Sample with an STL container.
Definition: itkListSample.h:52
unsigned int MeasurementVectorSizeType
Definition: itkSample.h:94
Converts pixel data into a list of pairs in order to compute a cooccurrence Histogram.
DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override
Make a DataObject of the correct type to used as the specified output.
void PrintSelf(std::ostream &os, Indent indent) const override
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
class ITK_FORWARD_EXPORT ProcessObject
Definition: itkDataObject.h:41