ITK
4.1.0
Insight Segmentation and Registration Toolkit
|
00001 /*========================================================================= 00002 * 00003 * Copyright Insight Software Consortium 00004 * 00005 * Licensed under the Apache License, Version 2.0 (the "License"); 00006 * you may not use this file except in compliance with the License. 00007 * You may obtain a copy of the License at 00008 * 00009 * http://www.apache.org/licenses/LICENSE-2.0.txt 00010 * 00011 * Unless required by applicable law or agreed to in writing, software 00012 * distributed under the License is distributed on an "AS IS" BASIS, 00013 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. 00014 * See the License for the specific language governing permissions and 00015 * limitations under the License. 00016 * 00017 *=========================================================================*/ 00018 /*========================================================================= 00019 * 00020 * Portions of this file are subject to the VTK Toolkit Version 3 copyright. 00021 * 00022 * Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen 00023 * 00024 * For complete copyright, license and disclaimer of warranty information 00025 * please refer to the NOTICE file at the top of the ITK source tree. 00026 * 00027 *=========================================================================*/ 00028 #ifndef __itkImageRegion_h 00029 #define __itkImageRegion_h 00030 00031 #include "itkRegion.h" 00032 00033 #include "itkSize.h" 00034 #include "itkContinuousIndex.h" 00035 #include "vnl/vnl_math.h" 00036 00037 namespace itk 00038 { 00039 // Forward declaration of ImageBase so it can be declared a friend 00040 // (needed for PrintSelf mechanism) 00041 template< unsigned int VImageDimension > 00042 class ImageBase; 00043 00068 template< unsigned int VImageDimension > 00069 class ITK_EXPORT ImageRegion:public Region 00070 { 00071 public: 00073 typedef ImageRegion Self; 00074 typedef Region Superclass; 00075 00077 itkTypeMacro(ImageRegion, Region); 00078 00080 itkStaticConstMacro(ImageDimension, unsigned int, VImageDimension); 00081 00084 itkStaticConstMacro( SliceDimension, unsigned int, 00085 ( ImageDimension - ( ImageDimension > 1 ) ) ); 00086 00088 static unsigned int GetImageDimension() 00089 { return ImageDimension; } 00090 00092 typedef Index< itkGetStaticConstMacro(ImageDimension) > IndexType; 00093 typedef typename IndexType::IndexValueType IndexValueType; 00094 typedef IndexValueType IndexValueArrayType[ImageDimension]; 00095 typedef typename IndexType::OffsetType OffsetType; 00096 typedef typename OffsetType::OffsetValueType OffsetValueType; 00097 00099 typedef Size< itkGetStaticConstMacro(ImageDimension) > SizeType; 00100 typedef typename SizeType::SizeValueType SizeValueType; 00101 00103 typedef ImageRegion< itkGetStaticConstMacro(SliceDimension) > SliceRegion; 00104 00106 virtual typename Superclass::RegionType GetRegionType() const 00107 { return Superclass::ITK_STRUCTURED_REGION; } 00108 00111 ImageRegion(); 00112 00115 virtual ~ImageRegion(); 00116 00119 ImageRegion(const Self & region):Region(region), m_Index(region.m_Index), m_Size(region.m_Size) {} 00120 00123 ImageRegion(const IndexType & index, const SizeType & size) 00124 { m_Index = index; m_Size = size; } 00125 00129 ImageRegion(const SizeType & size) 00130 { m_Size = size; m_Index.Fill(0); } 00131 00134 void operator=(const Self & region) 00135 { m_Index = region.m_Index; m_Size = region.m_Size; } 00136 00138 void SetIndex(const IndexType & index) 00139 { m_Index = index; } 00140 00142 const IndexType & GetIndex() const 00143 { return m_Index; } 00144 00147 void SetSize(const SizeType & size) 00148 { m_Size = size; } 00149 00151 const SizeType & GetSize() const 00152 { return m_Size; } 00153 00156 void SetSize(unsigned long i, SizeValueType sze) 00157 { m_Size[i] = sze; } 00158 SizeValueType GetSize(unsigned long i) const 00159 { return m_Size[i]; } 00161 00164 void SetIndex(unsigned long i, IndexValueType sze) 00165 { m_Index[i] = sze; } 00166 IndexValueType GetIndex(unsigned long i) const 00167 { return m_Index[i]; } 00169 00171 IndexType GetUpperIndex() const; 00172 00174 void SetUpperIndex( const IndexType & idx ); 00175 00177 bool 00178 operator==(const Self & region) const 00179 { 00180 bool same = 1; 00181 00182 same = ( m_Index == region.m_Index ); 00183 same = same && ( m_Size == region.m_Size ); 00184 return same; 00185 } 00186 00188 bool 00189 operator!=(const Self & region) const 00190 { 00191 bool same = 1; 00192 00193 same = ( m_Index == region.m_Index ); 00194 same = same && ( m_Size == region.m_Size ); 00195 return !same; 00196 } 00197 00199 bool 00200 IsInside(const IndexType & index) const 00201 { 00202 for ( unsigned int i = 0; i < ImageDimension; i++ ) 00203 { 00204 if ( index[i] < m_Index[i] ) 00205 { 00206 return false; 00207 } 00208 if ( index[i] >= ( m_Index[i] + static_cast< IndexValueType >( m_Size[i] ) ) ) 00209 { 00210 return false; 00211 } 00212 } 00213 return true; 00214 } 00216 00221 template< typename TCoordRepType > 00222 bool 00223 IsInside(const ContinuousIndex< TCoordRepType, VImageDimension > & index) const 00224 { 00225 for ( unsigned int i = 0; i < ImageDimension; i++ ) 00226 { 00227 if ( Math::RoundHalfIntegerUp< IndexValueType >(index[i]) < static_cast< IndexValueType >( m_Index[i] ) ) 00228 { 00229 return false; 00230 } 00231 // bound is the last valid pixel location 00232 const TCoordRepType bound = static_cast< TCoordRepType >( 00233 m_Index[i] + m_Size[i] - 0.5 ); 00235 00236 /* Note for NaN: test using negation of a positive test in order 00237 * to always evaluate to true (and thus return false) when index[i] 00238 * is NaN. The cast above to integer via RoundHalfIntegerUp will cast 00239 * NaN into a platform-dependent value (large negative, -1 or large 00240 * positive, empirically). Thus this test here is relied on 00241 * to 'catch' NaN's. */ 00242 if ( ! (index[i] <= bound) ) 00243 { 00244 return false; 00245 } 00246 } 00247 return true; 00248 } 00249 00254 bool 00255 IsInside(const Self & region) const 00256 { 00257 IndexType beginCorner = region.GetIndex(); 00258 00259 if ( !this->IsInside(beginCorner) ) 00260 { 00261 return false; 00262 } 00263 IndexType endCorner; 00264 SizeType size = region.GetSize(); 00265 for ( unsigned int i = 0; i < ImageDimension; i++ ) 00266 { 00267 endCorner[i] = beginCorner[i] + static_cast< OffsetValueType >( size[i] ) - 1; 00268 } 00269 if ( !this->IsInside(endCorner) ) 00270 { 00271 return false; 00272 } 00273 return true; 00274 } 00275 00278 SizeValueType GetNumberOfPixels() const; 00279 00283 void PadByRadius(OffsetValueType radius); 00284 00285 void PadByRadius(const IndexValueArrayType radius); 00286 00287 void PadByRadius(const SizeType & radius); 00288 00293 bool ShrinkByRadius(OffsetValueType radius); 00294 00295 bool ShrinkByRadius(const IndexValueArrayType radius); 00296 00297 bool ShrinkByRadius(const SizeType & radius); 00298 00303 bool Crop(const Self & region); 00304 00308 SliceRegion Slice(const unsigned long dim) const; 00309 00310 protected: 00315 virtual void PrintSelf(std::ostream & os, Indent indent) const; 00316 00317 private: 00318 IndexType m_Index; 00319 SizeType m_Size; 00320 00322 friend class ImageBase< VImageDimension >; 00323 }; 00324 00325 template< unsigned int VImageDimension > 00326 std::ostream & operator<<(std::ostream & os, const ImageRegion< VImageDimension > & region); 00327 } // end namespace itk 00328 00329 // Define instantiation macro for this template. 00330 #define ITK_TEMPLATE_ImageRegion(_, EXPORT, TypeX, TypeY) \ 00331 namespace itk \ 00332 { \ 00333 _( 1 ( class EXPORT ImageRegion< ITK_TEMPLATE_1 TypeX > ) ) \ 00334 _( 1 ( EXPORT std::ostream & operator<<(std::ostream &, \ 00335 const ImageRegion< ITK_TEMPLATE_1 TypeX > &) ) ) \ 00336 namespace Templates \ 00337 { \ 00338 typedef ImageRegion< ITK_TEMPLATE_1 TypeX > ImageRegion##TypeY; \ 00339 } \ 00340 } 00341 00342 #if ITK_TEMPLATE_EXPLICIT 00343 //template <unsigned int VImageDimension> const unsigned int 00344 // itk::ImageRegion<VImageDimension>::ImageDimension; 00345 //.template <unsigned int VImageDimension> const unsigned int 00346 // itk::ImageRegion<VImageDimension>::SliceDimension; 00347 #include "Templates/itkImageRegion+-.h" 00348 #endif 00349 00350 #if ITK_TEMPLATE_TXX 00351 #include "itkImageRegion.hxx" 00352 #endif 00353 00354 #endif 00355