ITK  4.2.0
Insight Segmentation and Registration Toolkit
itkImageRegion.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 /*=========================================================================
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 "vnl/vnl_math.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 ImageBase;
43 
68 template< unsigned int VImageDimension >
69 class ITK_EXPORT ImageRegion:public Region
70 {
71 public:
73  typedef ImageRegion Self;
74  typedef Region Superclass;
75 
77  itkTypeMacro(ImageRegion, Region);
78 
80  itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension);
81 
84  itkStaticConstMacro( SliceDimension, unsigned int,
85  ( ImageDimension - ( ImageDimension > 1 ) ) );
86 
88  static unsigned int GetImageDimension()
89  { return ImageDimension; }
90 
94  typedef IndexValueType IndexValueArrayType[ImageDimension];
97  typedef OffsetValueType OffsetTableType[ImageDimension+1];
98 
102 
105 
107  virtual typename Superclass::RegionType GetRegionType() const
108  { return Superclass::ITK_STRUCTURED_REGION; }
109 
112  ImageRegion();
113 
116  virtual ~ImageRegion();
117 
120  ImageRegion(const Self & region):Region(region), m_Index(region.m_Index), m_Size(region.m_Size) {}
121 
124  ImageRegion(const IndexType & index, const SizeType & size)
125  { m_Index = index; m_Size = size; }
126 
130  ImageRegion(const SizeType & size)
131  { m_Size = size; m_Index.Fill(0); }
132 
135  void operator=(const Self & region)
136  { m_Index = region.m_Index; m_Size = region.m_Size; }
137 
139  void SetIndex(const IndexType & index)
140  { m_Index = index; }
141 
143  const IndexType & GetIndex() const
144  { return m_Index; }
145 
148  void SetSize(const SizeType & size)
149  { m_Size = size; }
150 
152  const SizeType & GetSize() const
153  { return m_Size; }
154 
157  void SetSize(unsigned long i, SizeValueType sze)
158  { m_Size[i] = sze; }
159  SizeValueType GetSize(unsigned long i) const
160  { return m_Size[i]; }
162 
165  void SetIndex(unsigned long i, IndexValueType sze)
166  { m_Index[i] = sze; }
167  IndexValueType GetIndex(unsigned long i) const
168  { return m_Index[i]; }
170 
172  IndexType GetUpperIndex() const;
173 
175  void SetUpperIndex( const IndexType & idx );
176 
178  void ComputeOffsetTable(OffsetTableType offsetTable) const;
179 
181  bool
182  operator==(const Self & region) const
183  {
184  bool same = 1;
185 
186  same = ( m_Index == region.m_Index );
187  same = same && ( m_Size == region.m_Size );
188  return same;
189  }
190 
192  bool
193  operator!=(const Self & region) const
194  {
195  bool same = 1;
196 
197  same = ( m_Index == region.m_Index );
198  same = same && ( m_Size == region.m_Size );
199  return !same;
200  }
201 
203  bool
204  IsInside(const IndexType & index) const
205  {
206  for ( unsigned int i = 0; i < ImageDimension; i++ )
207  {
208  if ( index[i] < m_Index[i] )
209  {
210  return false;
211  }
212  if ( index[i] >= ( m_Index[i] + static_cast< IndexValueType >( m_Size[i] ) ) )
213  {
214  return false;
215  }
216  }
217  return true;
218  }
220 
225  template< typename TCoordRepType >
226  bool
228  {
229  for ( unsigned int i = 0; i < ImageDimension; i++ )
230  {
231  if ( Math::RoundHalfIntegerUp< IndexValueType >(index[i]) < static_cast< IndexValueType >( m_Index[i] ) )
232  {
233  return false;
234  }
235  // bound is the last valid pixel location
236  const TCoordRepType bound = static_cast< TCoordRepType >(
237  m_Index[i] + m_Size[i] - 0.5 );
239 
240  /* Note for NaN: test using negation of a positive test in order
241  * to always evaluate to true (and thus return false) when index[i]
242  * is NaN. The cast above to integer via RoundHalfIntegerUp will cast
243  * NaN into a platform-dependent value (large negative, -1 or large
244  * positive, empirically). Thus this test here is relied on
245  * to 'catch' NaN's. */
246  if ( ! (index[i] <= bound) )
247  {
248  return false;
249  }
250  }
251  return true;
252  }
253 
258  bool
259  IsInside(const Self & region) const
260  {
261  IndexType beginCorner = region.GetIndex();
262 
263  if ( !this->IsInside(beginCorner) )
264  {
265  return false;
266  }
267  IndexType endCorner;
268  SizeType size = region.GetSize();
269  for ( unsigned int i = 0; i < ImageDimension; i++ )
270  {
271  endCorner[i] = beginCorner[i] + static_cast< OffsetValueType >( size[i] ) - 1;
272  }
273  if ( !this->IsInside(endCorner) )
274  {
275  return false;
276  }
277  return true;
278  }
279 
282  SizeValueType GetNumberOfPixels() const;
283 
287  void PadByRadius(OffsetValueType radius);
288 
289  void PadByRadius(const IndexValueArrayType radius);
290 
291  void PadByRadius(const SizeType & radius);
292 
297  bool ShrinkByRadius(OffsetValueType radius);
298 
299  bool ShrinkByRadius(const IndexValueArrayType radius);
300 
301  bool ShrinkByRadius(const SizeType & radius);
302 
307  bool Crop(const Self & region);
308 
312  SliceRegion Slice(const unsigned long dim) const;
313 
314 protected:
319  virtual void PrintSelf(std::ostream & os, Indent indent) const;
320 
321 private:
324 
326  friend class ImageBase< VImageDimension >;
327 };
328 
329 template< unsigned int VImageDimension >
330 std::ostream & operator<<(std::ostream & os, const ImageRegion< VImageDimension > & region);
331 } // end namespace itk
332 
333 // Define instantiation macro for this template.
334 #define ITK_TEMPLATE_ImageRegion(_, EXPORT, TypeX, TypeY) \
335  namespace itk \
336  { \
337  _( 1 ( class EXPORT ImageRegion< ITK_TEMPLATE_1 TypeX > ) ) \
338  _( 1 ( EXPORT std::ostream & operator<<(std::ostream &, \
339  const ImageRegion< ITK_TEMPLATE_1 TypeX > &) ) ) \
340  namespace Templates \
341  { \
342  typedef ImageRegion< ITK_TEMPLATE_1 TypeX > ImageRegion##TypeY; \
343  } \
344  }
345 
346 #if ITK_TEMPLATE_EXPLICIT
347 //template <unsigned int VImageDimension> const unsigned int
348 // itk::ImageRegion<VImageDimension>::ImageDimension;
349 //.template <unsigned int VImageDimension> const unsigned int
350 // itk::ImageRegion<VImageDimension>::SliceDimension;
351 #include "Templates/itkImageRegion+-.h"
352 #endif
353 
354 #if ITK_TEMPLATE_TXX
355 #include "itkImageRegion.hxx"
356 #endif
357 
358 #endif
359