ITK  6.0.0
Insight Toolkit
itkImageRegion.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 /*=========================================================================
19  *
20  * Portions of this file are subject to the VTK Toolkit Version 3 copyright.
21  *
22  * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
23  *
24  * For complete copyright, license and disclaimer of warranty information
25  * please refer to the NOTICE file at the top of the ITK source tree.
26  *
27  *=========================================================================*/
28 #ifndef itkImageRegion_h
29 #define itkImageRegion_h
30 
31 #include "itkRegion.h"
32 
33 #include "itkSize.h"
34 #include "itkContinuousIndex.h"
35 #include "itkMath.h"
36 #include <type_traits> // For conditional and integral_constant.
37 #include <utility> // For tuple_element and tuple_size.
38 
39 // Macro added to each `ImageRegion` member function that overrides a virtual member function of `Region`, when legacy
40 // support is enabled. Without legacy support, `ImageRegion` will no longer inherit from `Region`, so then those
41 // `ImageRegion` member functions will no longer override.
42 #ifdef ITK_LEGACY_REMOVE
43 # define itkRegionOverrideMacro // nothing
44 #else
45 # define itkRegionOverrideMacro override
46 #endif
47 
48 namespace itk
49 {
50 // Forward declaration of ImageBase so it can be declared a friend
51 // (needed for PrintSelf mechanism)
52 template <unsigned int VImageDimension>
53 class ITK_TEMPLATE_EXPORT ImageBase;
54 
79 template <unsigned int VImageDimension>
80 class ITK_TEMPLATE_EXPORT ImageRegion final
81 #ifndef ITK_LEGACY_REMOVE
82  // This inheritance is only there when legacy support is enabled.
83  : public Region
84 #endif
85 {
86 public:
88  using Self = ImageRegion;
89 
90 #ifndef ITK_LEGACY_REMOVE
91  using Superclass = Region;
92 #endif
93 
95  const char *
97  {
98  return "ImageRegion";
99  }
100 
102  static constexpr unsigned int ImageDimension = VImageDimension;
103 
106  static constexpr unsigned int SliceDimension = ImageDimension - (ImageDimension > 1);
107 
109  static constexpr unsigned int
111  {
112  return ImageDimension;
113  }
114 
120  using IndexValueArrayType = IndexValueType[ImageDimension];
121  using OffsetTableType = OffsetValueType[ImageDimension + 1];
122 
126 
129 
133  {
135  }
136 
138  void
139  Print(std::ostream & os, Indent indent = 0) const itkRegionOverrideMacro;
140 
144  ImageRegion() noexcept = default;
145 
148  ~ImageRegion() itkRegionOverrideMacro = default;
149 
152  ImageRegion(const Self &) noexcept = default;
153 
157  ImageRegion(const IndexType & index, const SizeType & size) noexcept
158  : // Note: Use parentheses instead of curly braces to initialize data members,
159  // to avoid AppleClang 6.0.0.6000056 compile errors, "no viable conversion..."
160  m_Index(index)
161  , m_Size(size)
162  {}
163 
168  ImageRegion(const SizeType & size) noexcept
169  : m_Size(size)
170  {
171  // Note: m_Index is initialized by its C++11 default member initializer.
172  }
173 
176  Self &
177  operator=(const Self &) noexcept = default;
178 
180  void
181  SetIndex(const IndexType & index)
182  {
183  m_Index = index;
184  }
185 
187  const IndexType &
188  GetIndex() const
189  {
190  return m_Index;
191  }
192  IndexType &
194  {
195  return m_Index;
196  }
201  void
202  SetSize(const SizeType & size)
203  {
204  m_Size = size;
205  }
206 
208  const SizeType &
209  GetSize() const
210  {
211  return m_Size;
212  }
213  SizeType &
215  {
216  return m_Size;
217  }
222  void
223  SetSize(unsigned int i, SizeValueType sze)
224  {
225  m_Size[i] = sze;
226  }
228  GetSize(unsigned int i) const
229  {
230  return m_Size[i];
231  }
236  void
237  SetIndex(unsigned int i, IndexValueType sze)
238  {
239  m_Index[i] = sze;
240  }
242  GetIndex(unsigned int i) const
243  {
244  return m_Index[i];
245  }
249  IndexType
250  GetUpperIndex() const;
251 
253  void
254  SetUpperIndex(const IndexType & idx);
255 
257  void
258  ComputeOffsetTable(OffsetTableType offsetTable) const;
259 
261  bool
262  operator==(const Self & region) const noexcept
263  {
264  return (m_Index == region.m_Index) && (m_Size == region.m_Size);
265  }
266 
267  ITK_UNEQUAL_OPERATOR_MEMBER_FUNCTION(Self);
268 
270  bool
271  IsInside(const IndexType & index) const
272  {
273  for (unsigned int i = 0; i < ImageDimension; ++i)
274  {
275  if (index[i] < m_Index[i] || index[i] >= m_Index[i] + static_cast<IndexValueType>(m_Size[i]))
276  {
277  return false;
278  }
279  }
280  return true;
281  }
288  template <typename TCoordinate>
289  bool
291  {
292  constexpr TCoordinate half = 0.5;
293  for (unsigned int i = 0; i < ImageDimension; ++i)
294  {
295  // Use negation of tests so that index[i]==NaN leads to returning false.
296  if (!(index[i] >= m_Index[i] - half && index[i] <= (m_Index[i] + static_cast<IndexValueType>(m_Size[i])) - half))
297  {
298  return false;
299  }
300  }
301  return true;
302  }
309  bool
310  IsInside(const Self & otherRegion) const
311  {
312  const auto & otherIndex = otherRegion.m_Index;
313  const auto & otherSize = otherRegion.m_Size;
316  for (unsigned int i = 0; i < ImageDimension; ++i)
317  {
318  if (otherIndex[i] < m_Index[i] || otherSize[i] == 0 ||
319  otherIndex[i] + static_cast<IndexValueType>(otherSize[i]) >
320  m_Index[i] + static_cast<IndexValueType>(m_Size[i]))
321  {
322  return false;
323  }
324  }
325  return true;
326  }
327 
331  GetNumberOfPixels() const;
332 
336  void
337  PadByRadius(OffsetValueType radius);
338 
339  void
340  PadByRadius(const IndexValueArrayType radius);
341 
342  void
343  PadByRadius(const SizeType & radius);
344 
349  bool
350  ShrinkByRadius(OffsetValueType radius);
351 
352  bool
353  ShrinkByRadius(const IndexValueArrayType radius);
354 
355  bool
356  ShrinkByRadius(const SizeType & radius);
357 
362  bool
363  Crop(const Self & region);
364 
368  SliceRegion
369  Slice(const unsigned int dim) const;
370 
373  template <size_t VTupleIndex>
374  [[nodiscard]] auto &
375  get()
376  {
377  if constexpr (VTupleIndex == 0)
378  {
379  return m_Index;
380  }
381  else
382  {
383  static_assert(VTupleIndex == 1);
384  return m_Size;
385  }
386  }
390  template <size_t VTupleIndex>
391  [[nodiscard]] const auto &
392  get() const
393  {
394  if constexpr (VTupleIndex == 0)
395  {
396  return m_Index;
397  }
398  else
399  {
400  static_assert(VTupleIndex == 1);
401  return m_Size;
402  }
403  }
406 protected:
411  void
412  PrintSelf(std::ostream & os, Indent indent) const itkRegionOverrideMacro;
413 
414 private:
415  IndexType m_Index = { { 0 } };
416  SizeType m_Size = { { 0 } };
417 
419  friend class ImageBase<VImageDimension>;
420 };
421 
422 
423 // Deduction guide to avoid compiler warnings (-wctad-maybe-unsupported) when using class template argument deduction.
424 template <unsigned int VImageDimension>
426 
427 
428 template <unsigned int VImageDimension>
429 std::ostream &
430 operator<<(std::ostream & os, const ImageRegion<VImageDimension> & region);
431 } // end namespace itk
432 
433 
434 namespace std
435 {
436 #if defined(__clang__) && defined(__apple_build_version__) && (__clang_major__ <= 10)
437 # pragma clang diagnostic push
438 // Old AppleClang 10.0.0 (Xcode 10.1, newest on macOS 10.13) produced some unimportant warnings, like:
439 // "warning: 'tuple_size' defined as a struct template here but previously declared as a class template"
440 # pragma clang diagnostic ignored "-Wmismatched-tags"
441 #endif
442 
443 // NOLINTBEGIN(cert-dcl58-cpp)
444 // Locally suppressed the following warning from Clang-Tidy (LLVM 17.0.1), as it appears undeserved.
445 // > warning: modification of 'std' namespace can result in undefined behavior [cert-dcl58-cpp]
446 
454 template <unsigned int VImageDimension>
455 struct tuple_size<itk::ImageRegion<VImageDimension>> : integral_constant<size_t, 2>
456 {};
457 
459 template <size_t VTupleIndex, unsigned int VImageDimension>
460 struct tuple_element<VTupleIndex, itk::ImageRegion<VImageDimension>>
461  : conditional<VTupleIndex == 0, itk::Index<VImageDimension>, itk::Size<VImageDimension>>
462 {
463  static_assert(VTupleIndex < tuple_size_v<itk::ImageRegion<VImageDimension>>);
464 };
465 
466 // NOLINTEND(cert-dcl58-cpp)
467 
468 #if defined(__clang__) && defined(__apple_build_version__) && (__clang_major__ <= 10)
469 # pragma clang diagnostic pop
470 #endif
471 } // namespace std
472 
473 #undef itkRegionOverrideMacro
474 
475 #ifndef ITK_MANUAL_INSTANTIATION
476 # include "itkImageRegion.hxx"
477 #endif
478 
479 #endif
itk::ImageRegion::GetModifiableIndex
IndexType & GetModifiableIndex()
Definition: itkImageRegion.h:193
itk::ImageRegion< VDimension >::OffsetType
typename IndexType::OffsetType OffsetType
Definition: itkImageRegion.h:118
itkRegionOverrideMacro
#define itkRegionOverrideMacro
Definition: itkImageRegion.h:43
itk::ImageRegion< VDimension >::IndexValueType
typename IndexType::IndexValueType IndexValueType
Definition: itkImageRegion.h:117
itk::Index< VImageDimension >
itk::ImageRegion< VDimension >::OffsetTableType
OffsetValueType[ImageDimension+1] OffsetTableType
Definition: itkImageRegion.h:121
itk::ImageRegion::IsInside
bool IsInside(const Self &otherRegion) const
Definition: itkImageRegion.h:310
itkContinuousIndex.h
itk::Size< VImageDimension >
itk::ImageRegion< VDimension >::IndexValueArrayType
IndexValueType[ImageDimension] IndexValueArrayType
Definition: itkImageRegion.h:120
itk::ImageRegion::m_Size
SizeType m_Size
Definition: itkImageRegion.h:416
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:114
itk::ImageRegion::SetSize
void SetSize(unsigned int i, SizeValueType sze)
Definition: itkImageRegion.h:223
itk::ImageRegion::GetIndex
const IndexType & GetIndex() const
Definition: itkImageRegion.h:188
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:80
itkRegion.h
itk::ObjectEnums::RegionEnum::ITK_STRUCTURED_REGION
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::ImageRegion::GetSize
const SizeType & GetSize() const
Definition: itkImageRegion.h:209
itk::operator<<
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)
itk::ImageBase
class ITK_TEMPLATE_EXPORT ImageBase
Definition: itkImageHelper.h:50
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageRegion
ImageRegion(const Index< VImageDimension > &, const Size< VImageDimension > &) -> ImageRegion< VImageDimension >
itk::ImageRegion::SetIndex
void SetIndex(unsigned int i, IndexValueType sze)
Definition: itkImageRegion.h:237
itk::IndexValueType
long IndexValueType
Definition: itkIntTypes.h:93
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageRegion::operator==
bool operator==(const Self &region) const noexcept
Definition: itkImageRegion.h:262
itk::ImageRegion::GetSize
SizeValueType GetSize(unsigned int i) const
Definition: itkImageRegion.h:228
itk::ImageRegion::IsInside
bool IsInside(const IndexType &index) const
Definition: itkImageRegion.h:271
itk::ImageRegion::m_Index
IndexType m_Index
Definition: itkImageRegion.h:415
itk::ImageRegion::GetRegionType
Region::RegionEnum GetRegionType() const itkRegionOverrideMacro
Definition: itkImageRegion.h:132
itk::ImageRegion::get
auto & get()
Definition: itkImageRegion.h:375
itk::Region
A region represents some portion or piece of data.
Definition: itkRegion.h:65
itk::ImageRegion::GetImageDimension
static constexpr unsigned int GetImageDimension()
Definition: itkImageRegion.h:110
itk::OffsetValueType
long OffsetValueType
Definition: itkIntTypes.h:97
itk::ImageRegion::GetNameOfClass
const char * GetNameOfClass() const itkRegionOverrideMacro
Definition: itkImageRegion.h:96
itk::ObjectEnums::RegionEnum
RegionEnum
Definition: itkCommonEnums.h:253
itk::Offset
Represent a n-dimensional offset between two n-dimensional indexes of n-dimensional image.
Definition: itkOffset.h:66
itk::ImageRegion::IsInside
bool IsInside(const ContinuousIndex< TCoordinate, VImageDimension > &index) const
Definition: itkImageRegion.h:290
itk::ImageRegion::get
const auto & get() const
Definition: itkImageRegion.h:392
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnatomicalOrientation.h:29
itk::ContinuousIndex
A templated class holding a point in n-Dimensional image space.
Definition: itkContinuousIndex.h:46
itk::ImageRegion::SetIndex
void SetIndex(const IndexType &index)
Definition: itkImageRegion.h:181
itk::ImageRegion::GetIndex
IndexValueType GetIndex(unsigned int i) const
Definition: itkImageRegion.h:242
itk::ImageRegion< VDimension >::OffsetValueType
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageRegion.h:119
itk::ImageRegion< VDimension >::SizeValueType
typename SizeType::SizeValueType SizeValueType
Definition: itkImageRegion.h:125
AddImageFilter
Definition: itkAddImageFilter.h:81
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:202
Superclass
BinaryGeneratorImageFilter< TInputImage1, TInputImage2, TOutputImage > Superclass
Definition: itkAddImageFilter.h:90
itkMath.h
itk::ImageRegion::GetModifiableSize
SizeType & GetModifiableSize()
Definition: itkImageRegion.h:214
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:86
itk::ImageRegion::ImageRegion
ImageRegion(const SizeType &size) noexcept
Definition: itkImageRegion.h:168
itkSize.h