ITK  5.4.0
Insight Toolkit
itkLabelGeometryImageFilter.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 itkLabelGeometryImageFilter_h
19#define itkLabelGeometryImageFilter_h
20
22#include "itkNumericTraits.h"
23#include "itkArray.h"
25#include <map>
26#include <vector>
27#include "vnl/algo/vnl_symmetric_eigensystem.h"
28#include "vnl/vnl_det.h"
29#include "itkMath.h"
30
31namespace itk
32{
80template <typename TLabelImage, typename TIntensityImage = TLabelImage>
81class ITK_TEMPLATE_EXPORT [[deprecated(
82 "This class contains known computational bugs. See class documentation for details.")]] LabelGeometryImageFilter
83 : public ImageToImageFilter<TLabelImage, TIntensityImage>
84{
85public:
86 ITK_DISALLOW_COPY_AND_MOVE(LabelGeometryImageFilter);
87
93
95 itkNewMacro(Self);
96
98 itkOverrideGetNameOfClassMacro(LabelGeometryImageFilter);
99
101 using IntensityImageType = TIntensityImage;
106 using PixelType = typename TIntensityImage::PixelType;
107
109 using LabelImageType = TLabelImage;
114 using LabelPixelType = typename TLabelImage::PixelType;
116
118 static constexpr unsigned int ImageDimension = TLabelImage::ImageDimension;
119
122
125
128
132
133 // using BoundingBoxVerticesType = itk::FixedArray<
134 // LabelPointType,std::pow(2.0,Self::ImageDimension)>;
135 using BoundingBoxVerticesType = std::vector<LabelPointType>;
136
139
142
144 using LabelsType = std::vector<LabelPixelType>;
145
147 using LabelIndicesType = std::vector<LabelIndexType>;
148
150 using VectorType = std::vector<double>;
151
153 using MatrixType = vnl_matrix<double>;
154
160 {
161 public:
162 // default constructor
164 {
165 // initialized to the default values
166 this->m_Label = 0;
167 this->m_Sum = RealType{};
168
169 const unsigned int imageDimension = Self::ImageDimension;
170
171 // m_BoundingBox.resize(imageDimension*2);
172 for (unsigned int i = 0; i < imageDimension * 2; i += 2)
173 {
176 }
177
178 m_BoundingBoxVolume = 0;
179 m_BoundingBoxSize.Fill(0);
180 m_PixelIndices.clear();
181 m_Centroid.Fill(0);
182 m_WeightedCentroid.Fill(0);
183 m_ZeroOrderMoment = 0;
184 m_FirstOrderRawMoments.Fill(0);
185 m_FirstOrderWeightedRawMoments.Fill(0);
186 m_Eigenvalues.resize(ImageDimension);
187 m_Eigenvalues.clear();
188 m_Eigenvectors.set_size(ImageDimension, ImageDimension);
189 m_Eigenvectors.fill(0);
190 m_AxesLength.Fill(0);
191 m_Eccentricity = 1;
192 m_Elongation = 1;
193 m_Orientation = 0;
194 LabelPointType emptyPoint;
195 emptyPoint.Fill(0);
196 unsigned int numberOfVertices = 1 << ImageDimension;
197 m_OrientedBoundingBoxVertices.resize(numberOfVertices, emptyPoint);
198 m_OrientedBoundingBoxVolume = 0;
199 m_OrientedBoundingBoxSize.Fill(0);
200 m_OrientedLabelImage = LabelImageType::New();
201 m_OrientedIntensityImage = IntensityImageType::New();
202 m_OrientedBoundingBoxOrigin.Fill(0);
203 m_RotationMatrix.set_size(ImageDimension, ImageDimension);
204 m_RotationMatrix.fill(0.0);
205
206 m_SecondOrderRawMoments.set_size(ImageDimension, ImageDimension);
207 m_SecondOrderCentralMoments.set_size(ImageDimension, ImageDimension);
208 for (unsigned int i = 0; i < ImageDimension; ++i)
209 {
210 for (unsigned int j = 0; j < ImageDimension; ++j)
211 {
212 m_SecondOrderRawMoments(i, j) = 0;
213 m_SecondOrderCentralMoments(i, j) = 0;
214 }
215 }
216 }
217
246 };
247
249 // Map from the label to the class storing all of the geometry information.
250 using MapType = std::map<LabelPixelType, LabelGeometry>;
251 using MapIterator = typename std::map<LabelPixelType, LabelGeometry>::iterator;
252 using MapConstIterator = typename std::map<LabelPixelType, LabelGeometry>::const_iterator;
253
254 // Macros for enabling the calculation of additional features.
255 itkGetMacro(CalculatePixelIndices, bool);
256 itkBooleanMacro(CalculatePixelIndices);
257 void SetCalculatePixelIndices(const bool value)
258 {
259 // CalculateOrientedBoundingBox, CalculateOrientedLabelImage, and
260 // CalculateOrientedIntensityImage all need CalculatePixelIndices to be
261 // turned
262 // on if they are turned on. So, CalculatePixelIndices cannot be
263 // turned off if any of these flags are turned on.
264 if (value == false)
265 {
266 if ((this->m_CalculateOrientedBoundingBox) || (this->m_CalculateOrientedLabelRegions) ||
267 (this->m_CalculateOrientedIntensityRegions))
268 {
269 // We cannot change the value, so return.
270 return;
271 }
272 }
273
274 if (this->m_CalculatePixelIndices != value)
275 {
276 this->m_CalculatePixelIndices = value;
277 this->Modified();
278 }
279 }
280
281 itkGetMacro(CalculateOrientedBoundingBox, bool);
282 itkBooleanMacro(CalculateOrientedBoundingBox);
283 void SetCalculateOrientedBoundingBox(const bool value)
284 {
285 if (this->m_CalculateOrientedBoundingBox != value)
286 {
287 this->m_CalculateOrientedBoundingBox = value;
288 this->Modified();
289 }
290
291 // CalculateOrientedBoundingBox needs
292 // CalculatePixelIndices to be turned on.
293 if (value)
294 {
295 this->SetCalculatePixelIndices(true);
296 }
297 }
298
299 itkGetMacro(CalculateOrientedLabelRegions, bool);
300 itkBooleanMacro(CalculateOrientedLabelRegions);
301 void SetCalculateOrientedLabelRegions(const bool value)
302 {
303 if (this->m_CalculateOrientedLabelRegions != value)
304 {
305 this->m_CalculateOrientedLabelRegions = value;
306 this->Modified();
307
308 // CalculateOrientedLabelImage needs
309 // CalculateOrientedBoundingBox to be turned on.
310 if (value)
311 {
312 SetCalculateOrientedBoundingBox(true);
313 }
314 }
315 }
316
317 itkGetMacro(CalculateOrientedIntensityRegions, bool);
318 itkBooleanMacro(CalculateOrientedIntensityRegions);
320 {
321 if (this->m_CalculateOrientedIntensityRegions != value)
322 {
323 this->m_CalculateOrientedIntensityRegions = value;
324 this->Modified();
325
326 // CalculateOrientedIntensityImage needs
327 // CalculateOrientedBoundingBox to be turned on.
328 if (value)
329 {
330 this->SetCalculateOrientedBoundingBox(true);
331 }
332 }
333 }
334
336 void SetIntensityInput(const TIntensityImage * input)
337 {
338 // Process object is not const-correct so the const casting is required.
339 this->SetNthInput(1, const_cast<TIntensityImage *>(input));
340 }
341
343 const TIntensityImage * GetIntensityInput() const
344 {
345 return static_cast<TIntensityImage *>(const_cast<DataObject *>(this->ProcessObject::GetInput(1)));
346 }
347
350 bool HasLabel(LabelPixelType label) const { return m_LabelGeometryMapper.find(label) != m_LabelGeometryMapper.end(); }
351
353 SizeValueType GetNumberOfObjects() const { return m_LabelGeometryMapper.size(); }
354
355 SizeValueType GetNumberOfLabels() const { return this->GetNumberOfObjects(); }
356
358 std::vector<LabelPixelType> GetLabels() const { return m_AllLabels; }
359
362
366
369
372
375
378
381
384
388
392
395
399
402
407
410
413
421
424
427
430
434
437
439 TLabelImage * GetOrientedLabelImage(LabelPixelType label) const;
440
443 TIntensityImage * GetOrientedIntensityImage(LabelPixelType label) const;
444
445#ifdef ITK_USE_CONCEPT_CHECKING
446 // Begin concept checking
447 itkConceptMacro(InputHasNumericTraitsCheck, (Concept::HasNumericTraits<PixelType>));
448 // End concept checking
449#endif
450
451protected:
453 ~LabelGeometryImageFilter() override = default;
454 void PrintSelf(std::ostream & os, Indent indent) const override;
455
456 void GenerateData() override;
457
458private:
459 bool CalculateOrientedBoundingBoxVertices(vnl_symmetric_eigensystem<double> eig, LabelGeometry & m_LabelGeometry);
460
461 bool m_CalculatePixelIndices{};
462 bool m_CalculateOrientedBoundingBox{};
463 bool m_CalculateOrientedLabelRegions{};
464 bool m_CalculateOrientedIntensityRegions{};
465
466 MapType m_LabelGeometryMapper{};
467 LabelGeometry m_LabelGeometry{};
468 LabelsType m_AllLabels{};
469}; // end of class
470
471} // end namespace itk
472
473#ifndef ITK_MANUAL_INSTANTIATION
474# include "itkLabelGeometryImageFilter.hxx"
475#endif
476
477#endif
Base class for all data objects in ITK.
SmartPointer< Self > Pointer
Base class for filters that take an image as input and produce an image as output.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
FixedArray< float, Self::ImageDimension > m_AxesLength
Given a label map and an optional intensity image, compute geometric features.
typename TLabelImage::RegionType LabelRegionType
void SetCalculateOrientedIntensityRegions(const bool value)
RealType GetMajorAxisLength(LabelPixelType label) const
const TIntensityImage * GetIntensityInput() const
RealType GetElongation(LabelPixelType label) const
typename std::map< LabelPixelType, LabelGeometry >::iterator MapIterator
LabelSizeType GetBoundingBoxSize(LabelPixelType label) const
VectorType GetEigenvalues(LabelPixelType label) const
typename TLabelImage::SizeType LabelSizeType
typename TIntensityImage::SizeType SizeType
typename TIntensityImage::Pointer InputImagePointer
RealType GetEccentricity(LabelPixelType label) const
void SetCalculateOrientedBoundingBox(const bool value)
typename TIntensityImage::PixelType PixelType
typename TLabelImage::Pointer LabelImagePointer
LabelPointType GetCentroid(LabelPixelType label) const
void PrintSelf(std::ostream &os, Indent indent) const override
BoundingBoxVerticesType GetOrientedBoundingBoxVertices(LabelPixelType label) const
std::vector< LabelPixelType > LabelsType
LabelIndicesType GetPixelIndices(LabelPixelType label) const
typename TLabelImage::PixelType LabelPixelType
bool CalculateOrientedBoundingBoxVertices(vnl_symmetric_eigensystem< double > eig, LabelGeometry &m_LabelGeometry)
std::vector< LabelPixelType > GetLabels() const
std::vector< LabelPointType > BoundingBoxVerticesType
RealType GetIntegratedIntensity(LabelPixelType label) const
void SetCalculatePixelIndices(const bool value)
SizeValueType GetVolume(LabelPixelType label) const
MatrixType GetRotationMatrix(LabelPixelType label) const
AxesLengthType GetAxesLength(LabelPixelType label) const
typename TIntensityImage::IndexType IndexType
RealType GetBoundingBoxVolume(LabelPixelType label) const
TLabelImage * GetOrientedLabelImage(LabelPixelType label) const
bool HasLabel(LabelPixelType label) const
LabelPointType GetWeightedCentroid(LabelPixelType label) const
TIntensityImage * GetOrientedIntensityImage(LabelPixelType label) const
typename TIntensityImage::RegionType RegionType
RealType GetOrientedBoundingBoxVolume(LabelPixelType label) const
RealType GetMinorAxisLength(LabelPixelType label) const
MatrixType GetEigenvectors(LabelPixelType label) const
BoundingBoxType GetBoundingBox(LabelPixelType label) const
typename TLabelImage::PointType LabelPointType
LabelPointType GetOrientedBoundingBoxSize(LabelPixelType label) const
std::vector< LabelIndexType > LabelIndicesType
typename TLabelImage::IndexType LabelIndexType
typename std::map< LabelPixelType, LabelGeometry >::const_iterator MapConstIterator
std::map< LabelPixelType, LabelGeometry > MapType
void SetIntensityInput(const TIntensityImage *input)
~LabelGeometryImageFilter() override=default
LabelPointType GetOrientedBoundingBoxOrigin(LabelPixelType label) const
void SetCalculateOrientedLabelRegions(const bool value)
RealType GetOrientation(LabelPixelType label) const
typename NumericTraits< PixelType >::RealType RealType
RegionType GetRegion(LabelPixelType label) const
Light weight base class for most itk classes.
Define additional traits for native types such as int or float.
static constexpr T NonpositiveMin()
static constexpr T max(const T &)
DataObject * GetInput(const DataObjectIdentifierType &key)
Return an input.
Decorates any "simple" data type (data types without smart pointers) with a DataObject API.
static Pointer New()
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:83