ITK  5.4.0
Insight Toolkit
itkShapeLabelObject.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 itkShapeLabelObject_h
19#define itkShapeLabelObject_h
20
21#include "itkLabelObject.h"
22#include "itkLabelMap.h"
23#include "itkMath.h"
24#include "itkAffineTransform.h"
25
26namespace itk
27{
42template <typename TLabel, unsigned int VImageDimension>
43class ITK_TEMPLATE_EXPORT ShapeLabelObject : public LabelObject<TLabel, VImageDimension>
44{
45public:
46 ITK_DISALLOW_COPY_AND_MOVE(ShapeLabelObject);
47
51 using LabelObjectType = typename Superclass::LabelObjectType;
55
57 itkNewMacro(Self);
58
60 itkOverrideGetNameOfClassMacro(ShapeLabelObject);
61
63
64 static constexpr unsigned int ImageDimension = VImageDimension;
65
66 using typename Superclass::IndexType;
67
68 using LabelType = TLabel;
69
70 using typename Superclass::LineType;
71
72 using typename Superclass::LengthType;
73
74 using typename Superclass::AttributeType;
75
77 static constexpr AttributeType NUMBER_OF_PIXELS = 100;
78
82 static constexpr AttributeType PHYSICAL_SIZE = 101;
83
87 static constexpr AttributeType CENTROID = 104;
88
89 static constexpr AttributeType BOUNDING_BOX = 105;
90
97 static constexpr AttributeType NUMBER_OF_PIXELS_ON_BORDER = 106;
98
106 static constexpr AttributeType PERIMETER_ON_BORDER = 107;
107
111 static constexpr AttributeType FERET_DIAMETER = 108;
112
114 static constexpr AttributeType PRINCIPAL_MOMENTS = 109;
115
117 static constexpr AttributeType PRINCIPAL_AXES = 110;
118
122 static constexpr AttributeType ELONGATION = 111;
123
125 static constexpr AttributeType PERIMETER = 112;
126
127 static constexpr AttributeType ROUNDNESS = 113;
128
132 static constexpr AttributeType EQUIVALENT_SPHERICAL_RADIUS = 114;
133
137 static constexpr AttributeType EQUIVALENT_SPHERICAL_PERIMETER = 115;
138
142 static constexpr AttributeType EQUIVALENT_ELLIPSOID_DIAMETER = 116;
143
144 static constexpr AttributeType FLATNESS = 117;
145
146 static constexpr AttributeType PERIMETER_ON_BORDER_RATIO = 118;
147
148
151 static constexpr AttributeType ORIENTED_BOUNDING_BOX_ORIGIN = 119;
152
153
161 static constexpr AttributeType ORIENTED_BOUNDING_BOX_SIZE = 120;
162
163 static AttributeType
164 GetAttributeFromName(const std::string & s)
165 {
166 if (s == "NumberOfPixels")
167 {
168 return NUMBER_OF_PIXELS;
169 }
170 else if (s == "PhysicalSize")
171 {
172 return PHYSICAL_SIZE;
173 }
174 else if (s == "Centroid")
175 {
176 return CENTROID;
177 }
178 else if (s == "BoundingBox")
179 {
180 return BOUNDING_BOX;
181 }
182 else if (s == "NumberOfPixelsOnBorder")
183 {
184 return NUMBER_OF_PIXELS_ON_BORDER;
185 }
186 else if (s == "PerimeterOnBorder")
187 {
188 return PERIMETER_ON_BORDER;
189 }
190 else if (s == "FeretDiameter")
191 {
192 return FERET_DIAMETER;
193 }
194 else if (s == "PrincipalMoments")
195 {
196 return PRINCIPAL_MOMENTS;
197 }
198 else if (s == "PrincipalAxes")
199 {
200 return PRINCIPAL_AXES;
201 }
202 else if (s == "Elongation")
203 {
204 return ELONGATION;
205 }
206 else if (s == "Perimeter")
207 {
208 return PERIMETER;
209 }
210 else if (s == "Roundness")
211 {
212 return ROUNDNESS;
213 }
214 else if (s == "EquivalentSphericalRadius")
215 {
216 return EQUIVALENT_SPHERICAL_RADIUS;
217 }
218 else if (s == "EquivalentSphericalPerimeter")
219 {
220 return EQUIVALENT_SPHERICAL_PERIMETER;
221 }
222 else if (s == "EquivalentEllipsoidDiameter")
223 {
224 return EQUIVALENT_ELLIPSOID_DIAMETER;
225 }
226 else if (s == "Flatness")
227 {
228 return FLATNESS;
229 }
230 else if (s == "PerimeterOnBorderRatio")
231 {
232 return PERIMETER_ON_BORDER_RATIO;
233 }
234 else if (s == "OrientedBoundingBoxOrigin")
235 {
236 return ORIENTED_BOUNDING_BOX_ORIGIN;
237 }
238 else if (s == "OrientedBoundingBoxSize")
239 {
240 return ORIENTED_BOUNDING_BOX_SIZE;
241 }
242 // can't recognize the name
243 return Superclass::GetAttributeFromName(s);
244 }
245
246 static std::string
248 {
249 std::string name;
250 switch (a)
251 {
252 case NUMBER_OF_PIXELS:
253 name = "NumberOfPixels";
254 break;
255 case PHYSICAL_SIZE:
256 name = "PhysicalSize";
257 break;
258 case CENTROID:
259 name = "Centroid";
260 break;
261 case BOUNDING_BOX:
262 name = "BoundingBox";
263 break;
264 case NUMBER_OF_PIXELS_ON_BORDER:
265 name = "NumberOfPixelsOnBorder";
266 break;
267 case PERIMETER_ON_BORDER:
268 name = "PerimeterOnBorder";
269 break;
270 case FERET_DIAMETER:
271 name = "FeretDiameter";
272 break;
273 case PRINCIPAL_MOMENTS:
274 name = "PrincipalMoments";
275 break;
276 case PRINCIPAL_AXES:
277 name = "PrincipalAxes";
278 break;
279 case ELONGATION:
280 name = "Elongation";
281 break;
282 case PERIMETER:
283 name = "Perimeter";
284 break;
285 case ROUNDNESS:
286 name = "Roundness";
287 break;
288 case EQUIVALENT_SPHERICAL_RADIUS:
289 name = "EquivalentSphericalRadius";
290 break;
291 case EQUIVALENT_SPHERICAL_PERIMETER:
292 name = "EquivalentSphericalPerimeter";
293 break;
294 case EQUIVALENT_ELLIPSOID_DIAMETER:
295 name = "EquivalentEllipsoidDiameter";
296 break;
297 case FLATNESS:
298 name = "Flatness";
299 break;
300 case PERIMETER_ON_BORDER_RATIO:
301 name = "PerimeterOnBorderRatio";
302 break;
303 case ORIENTED_BOUNDING_BOX_ORIGIN:
304 name = "OrientedBoundingBoxOrigin";
305 break;
306 case ORIENTED_BOUNDING_BOX_SIZE:
307 name = "OrientedBoundingBoxSize";
308 break;
309 default:
310 // can't recognize the name
311 name = Superclass::GetNameFromAttribute(a);
312 break;
313 }
314 return name;
315 }
316
318
320
322
324
325public:
327
329
331
334
335
336 const RegionType &
338 {
339 return m_BoundingBox;
340 }
341
342 void
344 {
345 m_BoundingBox = v;
346 }
347
348 const double &
350 {
351 return m_PhysicalSize;
352 }
353
354 void
355 SetPhysicalSize(const double v)
356 {
357 m_PhysicalSize = v;
358 }
359
360 const SizeValueType &
362 {
363 return m_NumberOfPixels;
364 }
365
366 void
368 {
369 m_NumberOfPixels = v;
370 }
371
372 const CentroidType &
374 {
375 return m_Centroid;
376 }
377
378 void
379 SetCentroid(const CentroidType & centroid)
380 {
381 m_Centroid = centroid;
382 }
383
384 const SizeValueType &
386 {
387 return m_NumberOfPixelsOnBorder;
388 }
389
390 void
392 {
393 m_NumberOfPixelsOnBorder = v;
394 }
395
396 const double &
398 {
399 return m_PerimeterOnBorder;
400 }
401
402 void
403 SetPerimeterOnBorder(const double v)
404 {
405 m_PerimeterOnBorder = v;
406 }
407
408 const double &
410 {
411 return m_FeretDiameter;
412 }
413
414 void
415 SetFeretDiameter(const double v)
416 {
417 m_FeretDiameter = v;
418 }
419
420 const VectorType &
422 {
423 return m_PrincipalMoments;
424 }
425
426 void
428 {
429 m_PrincipalMoments = v;
430 }
431
432 const MatrixType &
434 {
435 return m_PrincipalAxes;
436 }
437
438 void
440 {
441 m_PrincipalAxes = v;
442 }
443
444 const double &
446 {
447 return m_Elongation;
448 }
449
450 void
451 SetElongation(const double v)
452 {
453 m_Elongation = v;
454 }
455
456 const double &
458 {
459 return m_Perimeter;
460 }
461
462 void
463 SetPerimeter(const double v)
464 {
465 m_Perimeter = v;
466 }
467
468 const double &
470 {
471 return m_Roundness;
472 }
473
474 void
475 SetRoundness(const double v)
476 {
477 m_Roundness = v;
478 }
479
480 const double &
482 {
483 return m_EquivalentSphericalRadius;
484 }
485
486 void
488 {
489 m_EquivalentSphericalRadius = v;
490 }
491
492 const double &
494 {
495 return m_EquivalentSphericalPerimeter;
496 }
497
498 void
500 {
501 m_EquivalentSphericalPerimeter = v;
502 }
503
504 const VectorType &
506 {
507 return m_EquivalentEllipsoidDiameter;
508 }
509
510 void
512 {
513 m_EquivalentEllipsoidDiameter = v;
514 }
515
516 const double &
518 {
519 return m_Flatness;
520 }
521
522 void
523 SetFlatness(const double v)
524 {
525 m_Flatness = v;
526 }
527
528 const double &
530 {
531 return m_PerimeterOnBorderRatio;
532 }
533
534 void
536 {
537 m_PerimeterOnBorderRatio = v;
538 }
539
540 const OrientedBoundingBoxPointType &
542 {
543 return m_OrientedBoundingBoxOrigin;
544 }
545
546 void
548 {
549 m_OrientedBoundingBoxOrigin = v;
550 }
551
552 const OrientedBoundingBoxSizeType &
554 {
555 return m_OrientedBoundingBoxSize;
556 }
557
558 void
560 {
561 m_OrientedBoundingBoxSize = v;
562 }
563
564
565 // some helper methods - not really required, but really useful!
566
568 const RegionType &
569 GetRegion() const
570 {
571 return m_BoundingBox;
572 }
573
574
577 const OrientedBoundingBoxDirectionType &
579 {
580 return this->GetPrincipalAxes();
581 }
582
593 OrientedBoundingBoxVerticesType
595 {
596 const MatrixType obbToPhysical(this->GetOrientedBoundingBoxDirection().GetTranspose());
597
598
600
601 // Use binary index to map the vertices of the OBB to an array. For
602 // example, in 2D, binary counting will give[0,0], [0,1], [1,0],
603 // [1,1], which corresponds to [minX,minY], [minX,maxY],
604 // [maxX,minY], [maxX,maxY].
605 for (unsigned int i = 0; i < OrientedBoundingBoxVerticesType::Length; ++i)
606 {
607 constexpr unsigned int msb = 1 << (ImageDimension - 1);
609 for (unsigned int j = 0; j < ImageDimension; ++j)
610 {
611 if (i & msb >> j)
612 {
613 offset[j] = m_OrientedBoundingBoxSize[j];
614 }
615 else
616 {
617 offset[j] = 0;
618 }
619 }
620 vertices[i] = m_OrientedBoundingBoxOrigin + obbToPhysical * offset;
621 }
622 return vertices;
623 }
624
628
634 {
635 typename AffineTransformType::MatrixType matrix;
636 typename AffineTransformType::OffsetType offset;
637 for (unsigned int i = 0; i < VImageDimension; ++i)
638 {
639 offset[i] = m_Centroid[i];
640 for (unsigned int j = 0; j < VImageDimension; ++j)
641 {
642 matrix[j][i] = m_PrincipalAxes[i][j]; // Note the transposition
643 }
644 }
648
649 result->SetMatrix(matrix);
650 result->SetOffset(offset);
651
652 return result;
653 }
654
659 AffineTransformPointer
661 {
662 typename AffineTransformType::MatrixType matrix;
663 typename AffineTransformType::OffsetType offset;
664 for (unsigned int i = 0; i < VImageDimension; ++i)
665 {
666 offset[i] = m_Centroid[i];
667 for (unsigned int j = 0; j < VImageDimension; ++j)
668 {
669 matrix[j][i] = m_PrincipalAxes[i][j]; // Note the transposition
670 }
671 }
675 result->SetMatrix(matrix);
676 result->SetOffset(offset);
677
679 result->GetInverse(inverse);
680
681 return inverse;
682 }
683
684 template <typename TSourceLabelObject>
685 void
686 CopyAttributesFrom(const TSourceLabelObject * src)
687 {
688 Superclass::template CopyAttributesFrom<TSourceLabelObject>(src);
689
690 m_BoundingBox = src->GetBoundingBox();
691 m_NumberOfPixels = src->GetNumberOfPixels();
692 m_PhysicalSize = src->GetPhysicalSize();
693 m_Centroid = src->GetCentroid();
694 m_NumberOfPixelsOnBorder = src->GetNumberOfPixelsOnBorder();
695 m_PerimeterOnBorder = src->GetPerimeterOnBorder();
696 m_FeretDiameter = src->GetFeretDiameter();
697 m_PrincipalMoments = src->GetPrincipalMoments();
698 m_PrincipalAxes = src->GetPrincipalAxes();
699 m_Elongation = src->GetElongation();
700 m_Perimeter = src->GetPerimeter();
701 m_Roundness = src->GetRoundness();
702 m_EquivalentSphericalRadius = src->GetEquivalentSphericalRadius();
703 m_EquivalentSphericalPerimeter = src->GetEquivalentSphericalPerimeter();
704 m_EquivalentEllipsoidDiameter = src->GetEquivalentEllipsoidDiameter();
705 m_Flatness = src->GetFlatness();
706 m_PerimeterOnBorderRatio = src->GetPerimeterOnBorderRatio();
707 m_OrientedBoundingBoxOrigin = src->GetOrientedBoundingBoxOrigin();
708 m_OrientedBoundingBoxSize = src->GetOrientedBoundingBoxSize();
709 }
710
711 template <typename TSourceLabelObject>
712 void
713 CopyAllFrom(const TSourceLabelObject * src)
714 {
715 itkAssertOrThrowMacro((src != nullptr), "Null Pointer");
716 this->template CopyLinesFrom<TSourceLabelObject>(src);
717 this->template CopyAttributesFrom<TSourceLabelObject>(src);
718 }
719
720protected:
722 {
723 m_NumberOfPixels = 0;
724 m_PhysicalSize = 0;
725 m_Centroid.Fill(0);
726 m_NumberOfPixelsOnBorder = 0;
727 m_PerimeterOnBorder = 0;
728 m_FeretDiameter = 0;
729 m_PrincipalMoments.Fill(0);
730 m_PrincipalAxes.Fill(0);
731 m_Elongation = 0;
732 m_Perimeter = 0;
733 m_Roundness = 0;
734 m_EquivalentSphericalRadius = 0;
735 m_EquivalentSphericalPerimeter = 0;
736 m_EquivalentEllipsoidDiameter.Fill(0);
737 m_Flatness = 0;
738 m_PerimeterOnBorderRatio = 0;
739 m_OrientedBoundingBoxSize.Fill(0);
740 m_OrientedBoundingBoxOrigin.Fill(0);
741 }
742
743 void
744 PrintSelf(std::ostream & os, Indent indent) const override
745 {
746 Superclass::PrintSelf(os, indent);
747
748 os << indent << "BoundingBox: " << m_BoundingBox << std::endl;
749 os << indent
750 << "NumberOfPixels: " << static_cast<typename NumericTraits<SizeValueType>::PrintType>(m_NumberOfPixels)
751 << std::endl;
752 os << indent << "PhysicalSize: " << m_PhysicalSize << std::endl;
753 os << indent << "Centroid: " << static_cast<typename NumericTraits<CentroidType>::PrintType>(m_Centroid)
754 << std::endl;
755 os << indent << "NumberOfPixelsOnBorder: "
756 << static_cast<typename NumericTraits<SizeValueType>::PrintType>(m_NumberOfPixelsOnBorder) << std::endl;
757 os << indent << "PerimeterOnBorder: " << m_PerimeterOnBorder << std::endl;
758 os << indent << "FeretDiameter: " << m_FeretDiameter << std::endl;
759 os << indent << "PrincipalMoments: " << m_PrincipalMoments << std::endl;
760 os << indent << "PrincipalAxes: " << std::endl << m_PrincipalAxes;
761 os << indent << "Elongation: " << m_Elongation << std::endl;
762 os << indent << "Perimeter: " << m_Perimeter << std::endl;
763 os << indent << "Roundness: " << m_Roundness << std::endl;
764 os << indent << "EquivalentSphericalRadius: " << m_EquivalentSphericalRadius << std::endl;
765 os << indent << "EquivalentSphericalPerimeter: " << m_EquivalentSphericalPerimeter << std::endl;
766 os << indent << "EquivalentEllipsoidDiameter: " << m_EquivalentEllipsoidDiameter << std::endl;
767 os << indent << "Flatness: " << m_Flatness << std::endl;
768 os << indent << "PerimeterOnBorderRatio: " << m_PerimeterOnBorderRatio << std::endl;
769 os << indent << "OrientedBoundingBoxSize: "
770 << static_cast<typename NumericTraits<OrientedBoundingBoxSizeType>::PrintType>(m_OrientedBoundingBoxSize)
771 << std::endl;
772 os << indent << "OrientedBoundingBoxOrigin: "
773 << static_cast<typename NumericTraits<OrientedBoundingBoxPointType>::PrintType>(m_OrientedBoundingBoxOrigin)
774 << std::endl;
775 }
776
777private:
778 RegionType m_BoundingBox{};
779 SizeValueType m_NumberOfPixels{};
780 double m_PhysicalSize{};
781 CentroidType m_Centroid{};
782 SizeValueType m_NumberOfPixelsOnBorder{};
783 double m_PerimeterOnBorder{};
784 double m_FeretDiameter{};
785 VectorType m_PrincipalMoments{};
786 MatrixType m_PrincipalAxes{};
787 double m_Elongation{};
788 double m_Perimeter{};
789 double m_Roundness{};
790 double m_EquivalentSphericalRadius{};
791 double m_EquivalentSphericalPerimeter{};
792 VectorType m_EquivalentEllipsoidDiameter{};
793 double m_Flatness{};
794 double m_PerimeterOnBorderRatio{};
795
796 OrientedBoundingBoxSizeType m_OrientedBoundingBoxSize{};
797 OrientedBoundingBoxPointType m_OrientedBoundingBoxOrigin{};
798};
799} // end namespace itk
800
801#endif
Simulate a standard C array with copy semantics.
Definition: itkFixedArray.h:54
An image region represents a structured region of data.
Control indentation during Print() invocation.
Definition: itkIndent.h:50
Templated n-dimensional image to store labeled objects.
Definition: itkLabelMap.h:71
The base class for the representation of a labeled binary object in an image.
itk::SizeValueType SizeValueType
unsigned int AttributeType
Light weight base class for most itk classes.
A Label object to store the common attributes related to the shape of the object.
const double & GetEquivalentSphericalPerimeter() const
void SetEquivalentEllipsoidDiameter(const VectorType &v)
void CopyAllFrom(const TSourceLabelObject *src)
const VectorType & GetPrincipalMoments() const
void SetElongation(const double v)
const SizeValueType & GetNumberOfPixelsOnBorder() const
const double & GetPerimeter() const
const MatrixType & GetPrincipalAxes() const
void SetEquivalentSphericalRadius(const double v)
const RegionType & GetRegion() const
OrientedBoundingBoxVerticesType GetOrientedBoundingBoxVertices() const
const OrientedBoundingBoxPointType & GetOrientedBoundingBoxOrigin() const
void SetPrincipalMoments(const VectorType &v)
void SetPerimeterOnBorder(const double v)
const double & GetFlatness() const
void SetNumberOfPixels(const SizeValueType &v)
void SetNumberOfPixelsOnBorder(const SizeValueType &v)
const CentroidType & GetCentroid() const
static AttributeType GetAttributeFromName(const std::string &s)
const RegionType & GetBoundingBox() const
void SetPrincipalAxes(const MatrixType &v)
const double & GetEquivalentSphericalRadius() const
AffineTransformPointer GetPrincipalAxesToPhysicalAxesTransform() const
const double & GetFeretDiameter() const
void SetPerimeter(const double v)
void SetFeretDiameter(const double v)
const OrientedBoundingBoxSizeType & GetOrientedBoundingBoxSize() const
const OrientedBoundingBoxDirectionType & GetOrientedBoundingBoxDirection() const
const double & GetElongation() const
void SetOrientedBoundingBoxSize(const OrientedBoundingBoxSizeType &v)
void CopyAttributesFrom(const TSourceLabelObject *src)
const double & GetPhysicalSize() const
void SetPerimeterOnBorderRatio(const double v)
void SetBoundingBox(const RegionType &v)
const double & GetPerimeterOnBorder() const
const SizeValueType & GetNumberOfPixels() const
const VectorType & GetEquivalentEllipsoidDiameter() const
void PrintSelf(std::ostream &os, Indent indent) const override
const double & GetPerimeterOnBorderRatio() const
typename AffineTransformType::Pointer AffineTransformPointer
void SetEquivalentSphericalPerimeter(const double v)
const double & GetRoundness() const
void SetFlatness(const double v)
void SetCentroid(const CentroidType &centroid)
AffineTransformPointer GetPhysicalAxesToPrincipalAxesTransform() const
void SetOrientedBoundingBoxOrigin(const OrientedBoundingBoxPointType &v)
static std::string GetNameFromAttribute(const AttributeType &a)
void SetPhysicalSize(const double v)
void SetRoundness(const double v)
Implements a weak reference to an object.
static Pointer New()
SmartPointer< Self > Pointer
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:83