ITK  5.0.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:
98  ITK_DISALLOW_COPY_AND_ASSIGN(ScalarImageToCooccurrenceMatrixFilter);
99 
105 
108 
110  itkNewMacro(Self);
111 
112  using ImageType = TImageType;
113  using ImagePointer = typename ImageType::Pointer;
114  using ImageConstPointer = typename ImageType::ConstPointer;
115  using PixelType = typename ImageType::PixelType;
117  using RadiusType = typename ImageType::SizeType;
118  using OffsetType = typename ImageType::OffsetType;
122 
124 
129 
130  static constexpr unsigned int DefaultBinsPerAxis = 256;
131 
134  itkSetConstObjectMacro(Offsets, OffsetVector);
135  itkGetConstObjectMacro(Offsets, OffsetVector);
137 
138  void SetOffset(const OffsetType offset);
139 
141  itkSetMacro(NumberOfBinsPerAxis, unsigned int);
142  itkGetConstMacro(NumberOfBinsPerAxis, unsigned int);
144 
147  void SetPixelValueMinMax(PixelType min, PixelType max);
148 
149  itkGetConstMacro(Min, PixelType);
150  itkGetConstMacro(Max, PixelType);
151 
154  itkSetMacro(Normalize, bool);
155  itkGetConstMacro(Normalize, bool);
156  itkBooleanMacro(Normalize);
158 
160  using Superclass::SetInput;
161  void SetInput(const ImageType *image);
162 
163  const ImageType * GetInput() const;
164 
166  void SetMaskImage(const ImageType *image);
167 
168  const ImageType * GetMaskImage() const;
169 
171  const HistogramType * GetOutput() const;
172 
175  itkSetMacro(InsidePixelValue, PixelType);
176  itkGetConstMacro(InsidePixelValue, PixelType);
178 
179 protected:
181  ~ScalarImageToCooccurrenceMatrixFilter() override = default;
182  void PrintSelf(std::ostream & os, Indent indent) const override;
183 
184  virtual void FillHistogram(RadiusType radius, RegionType region);
185 
186  virtual void FillHistogramWithMask(RadiusType radius, RegionType region, const ImageType *maskImage);
187 
190 
192  using Superclass::MakeOutput;
193  DataObjectPointer MakeOutput(DataObjectPointerArraySizeType idx) override;
194 
196  void GenerateData() override;
197 
198 private:
199  void NormalizeHistogram();
200 
204 
205  unsigned int m_NumberOfBinsPerAxis;
209 
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.
DataObjectPointerArray::size_type DataObjectPointerArraySizeType
Define numeric traits for std::vector.
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...
class ITK_FORWARD_EXPORT ProcessObject
Definition: itkDataObject.h:40
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
typename Superclass::MeasurementVectorType MeasurementVectorType
Definition: itkHistogram.h:102
SmartPointer< Self > Pointer