ITK  5.2.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  * 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 /*=========================================================================
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 
37 namespace itk
38 {
39 // Forward declaration of ImageBase so it can be declared a friend
40 // (needed for PrintSelf mechanism)
41 template <unsigned int VImageDimension>
42 class ITK_TEMPLATE_EXPORT ImageBase;
43 
68 template <unsigned int VImageDimension>
69 class ITK_TEMPLATE_EXPORT ImageRegion final : public Region
70 {
71 public:
73  using Self = ImageRegion;
74  using Superclass = Region;
75 
77  itkTypeMacro(ImageRegion, Region);
78 
80  static constexpr unsigned int ImageDimension = VImageDimension;
81 
84  static constexpr unsigned int SliceDimension = ImageDimension - (ImageDimension > 1);
85 
87  static unsigned int
89  {
90  return ImageDimension;
91  }
92 
98  using IndexValueArrayType = IndexValueType[ImageDimension];
99  using OffsetTableType = OffsetValueType[ImageDimension + 1];
100 
104 
107 
110  GetRegionType() const override
111  {
112  return Superclass::RegionEnum::ITK_STRUCTURED_REGION;
113  }
114 
118  ImageRegion() ITK_NOEXCEPT = default;
119 
122  ~ImageRegion() override = default;
123 
126  ImageRegion(const Self &) ITK_NOEXCEPT = default;
127 
130  ImageRegion(const IndexType & index, const SizeType & size) ITK_NOEXCEPT
131  :
132  // Note: Use parentheses instead of curly braces to initialize data members,
133  // to avoid AppleClang 6.0.0.6000056 compile errors, "no viable conversion..."
134  m_Index(index)
135  , m_Size(size)
136  {}
137 
141  ImageRegion(const SizeType & size) ITK_NOEXCEPT : m_Size(size)
142  {
143  // Note: m_Index is initialized by its C++11 default member initializer.
144  }
145 
148  Self &
149  operator=(const Self &) ITK_NOEXCEPT = default;
150 
152  void
153  SetIndex(const IndexType & index)
154  {
155  m_Index = index;
156  }
157 
159  const IndexType &
160  GetIndex() const
161  {
162  return m_Index;
163  }
164  IndexType &
166  {
167  return m_Index;
168  }
170 
173  void
174  SetSize(const SizeType & size)
175  {
176  m_Size = size;
177  }
178 
180  const SizeType &
181  GetSize() const
182  {
183  return m_Size;
184  }
185  SizeType &
187  {
188  return m_Size;
189  }
191 
194  void
195  SetSize(unsigned int i, SizeValueType sze)
196  {
197  m_Size[i] = sze;
198  }
200  GetSize(unsigned int i) const
201  {
202  return m_Size[i];
203  }
205 
208  void
209  SetIndex(unsigned int i, IndexValueType sze)
210  {
211  m_Index[i] = sze;
212  }
214  GetIndex(unsigned int i) const
215  {
216  return m_Index[i];
217  }
219 
221  IndexType
222  GetUpperIndex() const;
223 
225  void
226  SetUpperIndex(const IndexType & idx);
227 
229  void
230  ComputeOffsetTable(OffsetTableType offsetTable) const;
231 
233  bool
234  operator==(const Self & region) const ITK_NOEXCEPT
235  {
236  return (m_Index == region.m_Index) && (m_Size == region.m_Size);
237  }
238 
240  bool
241  operator!=(const Self & region) const ITK_NOEXCEPT
242  {
243  return !(*this == region);
244  }
245 
247  bool
248  IsInside(const IndexType & index) const
249  {
250  for (unsigned int i = 0; i < ImageDimension; i++)
251  {
252  if (index[i] < m_Index[i])
253  {
254  return false;
255  }
256  if (index[i] >= (m_Index[i] + static_cast<IndexValueType>(m_Size[i])))
257  {
258  return false;
259  }
260  }
261  return true;
262  }
264 
269  template <typename TCoordRepType>
270  bool
272  {
273  for (unsigned int i = 0; i < ImageDimension; i++)
274  {
275  if (Math::RoundHalfIntegerUp<IndexValueType>(index[i]) < static_cast<IndexValueType>(m_Index[i]))
276  {
277  return false;
278  }
279  // bound is the last valid pixel location
280  const auto bound = static_cast<TCoordRepType>(m_Index[i] + m_Size[i] - 0.5);
282 
283  /* Note for NaN: test using negation of a positive test in order
284  * to always evaluate to true (and thus return false) when index[i]
285  * is NaN. The cast above to integer via RoundHalfIntegerUp will cast
286  * NaN into a platform-dependent value (large negative, -1 or large
287  * positive, empirically). Thus this test here is relied on
288  * to 'catch' NaN's. */
289  if (!(index[i] <= bound))
290  {
291  return false;
292  }
293  }
294  return true;
295  }
296 
301  bool
302  IsInside(const Self & region) const
303  {
304  IndexType beginCorner = region.GetIndex();
305 
306  if (!this->IsInside(beginCorner))
307  {
308  return false;
309  }
310  IndexType endCorner;
311  const SizeType & size = region.GetSize();
312  for (unsigned int i = 0; i < ImageDimension; i++)
313  {
314  endCorner[i] = beginCorner[i] + static_cast<OffsetValueType>(size[i]) - 1;
315  }
316  if (!this->IsInside(endCorner))
317  {
318  return false;
319  }
320  return true;
321  }
322 
326  GetNumberOfPixels() const;
327 
331  void
332  PadByRadius(OffsetValueType radius);
333 
334  void
335  PadByRadius(const IndexValueArrayType radius);
336 
337  void
338  PadByRadius(const SizeType & radius);
339 
344  bool
345  ShrinkByRadius(OffsetValueType radius);
346 
347  bool
348  ShrinkByRadius(const IndexValueArrayType radius);
349 
350  bool
351  ShrinkByRadius(const SizeType & radius);
352 
357  bool
358  Crop(const Self & region);
359 
363  SliceRegion
364  Slice(const unsigned int dim) const;
365 
366 protected:
371  void
372  PrintSelf(std::ostream & os, Indent indent) const override;
373 
374 private:
375  IndexType m_Index = { { 0 } };
376  SizeType m_Size = { { 0 } };
377 
379  friend class ImageBase<VImageDimension>;
380 };
381 
382 template <unsigned int VImageDimension>
383 std::ostream &
384 operator<<(std::ostream & os, const ImageRegion<VImageDimension> & region);
385 } // end namespace itk
386 
387 #ifndef ITK_MANUAL_INSTANTIATION
388 # include "itkImageRegion.hxx"
389 #endif
390 
391 #endif
itk::ImageRegion::GetModifiableIndex
IndexType & GetModifiableIndex()
Definition: itkImageRegion.h:165
itk::ImageRegion< VDimension >::OffsetType
typename IndexType::OffsetType OffsetType
Definition: itkImageRegion.h:96
itk::ImageRegion< VDimension >::IndexValueType
typename IndexType::IndexValueType IndexValueType
Definition: itkImageRegion.h:95
itk::Index< Self::ImageDimension >
itk::ImageRegion< VDimension >::OffsetTableType
OffsetValueType[ImageDimension+1] OffsetTableType
Definition: itkImageRegion.h:99
itkContinuousIndex.h
itk::Size< Self::ImageDimension >
itk::operator<<
std::ostream & operator<<(std::ostream &os, const Array< TValue > &arr)
Definition: itkArray.h:218
itk::ImageRegion< VDimension >::IndexValueArrayType
IndexValueType[ImageDimension] IndexValueArrayType
Definition: itkImageRegion.h:98
itk::ImageBase
Base class for templated image classes.
Definition: itkImageBase.h:105
itk::ImageRegion::SetSize
void SetSize(unsigned int i, SizeValueType sze)
Definition: itkImageRegion.h:195
itk::ImageRegion::GetIndex
const IndexType & GetIndex() const
Definition: itkImageRegion.h:160
itk::ImageRegion
An image region represents a structured region of data.
Definition: itkImageRegion.h:69
itkRegion.h
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::ImageRegion::GetSize
const SizeType & GetSize() const
Definition: itkImageRegion.h:181
itk::ImageBase
class ITK_TEMPLATE_EXPORT ImageBase
Definition: itkImageHelper.h:51
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::ImageRegion::SetIndex
void SetIndex(unsigned int i, IndexValueType sze)
Definition: itkImageRegion.h:209
itk::ObjectEnums::RegionEnum
RegionEnum
Definition: itkCommonEnums.h:258
itk::ImageRegion::IsInside
bool IsInside(const Self &region) const
Definition: itkImageRegion.h:302
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:234
itk::ImageRegion::GetSize
SizeValueType GetSize(unsigned int i) const
Definition: itkImageRegion.h:200
itk::ImageRegion::IsInside
bool IsInside(const IndexType &index) const
Definition: itkImageRegion.h:248
itk::ImageRegion::IsInside
bool IsInside(const ContinuousIndex< TCoordRepType, VImageDimension > &index) const
Definition: itkImageRegion.h:271
itk::Region
A region represents some portion or piece of data.
Definition: itkRegion.h:65
itk::ImageRegion::GetImageDimension
static unsigned int GetImageDimension()
Definition: itkImageRegion.h:88
itk::Offset
Represent a n-dimensional offset between two n-dimensional indexes of n-dimensional image.
Definition: itkOffset.h:67
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ContinuousIndex
A templated class holding a point in n-Dimensional image space.
Definition: itkContinuousIndex.h:46
itk::OffsetValueType
signed long OffsetValueType
Definition: itkIntTypes.h:94
itk::ImageRegion::GetRegionType
Superclass::RegionEnum GetRegionType() const override
Definition: itkImageRegion.h:110
itk::ImageRegion::SetIndex
void SetIndex(const IndexType &index)
Definition: itkImageRegion.h:153
itk::ImageRegion::GetIndex
IndexValueType GetIndex(unsigned int i) const
Definition: itkImageRegion.h:214
itk::IndexValueType
signed long IndexValueType
Definition: itkIntTypes.h:90
itk::ImageRegion< VDimension >::OffsetValueType
typename OffsetType::OffsetValueType OffsetValueType
Definition: itkImageRegion.h:97
itk::ImageRegion::operator!=
bool operator!=(const Self &region) const noexcept
Definition: itkImageRegion.h:241
itk::ImageRegion< VDimension >::SizeValueType
typename SizeType::SizeValueType SizeValueType
Definition: itkImageRegion.h:103
itk::ImageRegion::SetSize
void SetSize(const SizeType &size)
Definition: itkImageRegion.h:174
itkMath.h
itk::ImageRegion::GetModifiableSize
SizeType & GetModifiableSize()
Definition: itkImageRegion.h:186
itk::SizeValueType
unsigned long SizeValueType
Definition: itkIntTypes.h:83
itk::ImageRegion::ImageRegion
ImageRegion(const SizeType &size) noexcept
Definition: itkImageRegion.h:141
itkSize.h