Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

Iterators/ImageRegionIterator.cxx

See also:
ImageConstIterator

ConditionalConstIterator

ConstNeighborhoodIterator

ConstShapedNeighborhoodIterator

ConstSliceIterator

CorrespondenceDataStructureIterator

FloodFilledFunctionConditionalConstIterator

FloodFilledImageFunctionConditionalConstIterator

FloodFilledImageFunctionConditionalIterator

FloodFilledSpatialFunctionConditionalConstIterator

FloodFilledSpatialFunctionConditionalIterator

ImageConstIterator

ImageConstIteratorWithIndex

ImageIterator

ImageIteratorWithIndex

ImageLinearConstIteratorWithIndex

ImageLinearIteratorWithIndex

ImageRandomConstIteratorWithIndex

ImageRandomIteratorWithIndex

ImageRegionConstIterator

ImageRegionConstIteratorWithIndex

ImageRegionExclusionConstIteratorWithIndex

ImageRegionExclusionIteratorWithIndex

ImageRegionIterator

ImageRegionIteratorWithIndex

ImageRegionReverseConstIterator

ImageRegionReverseIterator

ImageReverseConstIterator

ImageReverseIterator

ImageSliceConstIteratorWithIndex

ImageSliceIteratorWithIndex

NeighborhoodIterator

PathConstIterator

PathIterator

ShapedNeighborhoodIterator

SliceIterator

ImageConstIteratorWithIndex

00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: ImageRegionIterator.cxx,v $
00005   Language:  C++
00006   Date:      $Date: 2005/02/08 03:59:00 $
00007   Version:   $Revision: 1.24 $
00008 
00009   Copyright (c) Insight Software Consortium. All rights reserved.
00010   See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
00011 
00012      This software is distributed WITHOUT ANY WARRANTY; without even 
00013      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
00014      PURPOSE.  See the above copyright notices for more information.
00015 
00016 =========================================================================*/
00017 #if defined(_MSC_VER)
00018 #pragma warning ( disable : 4786 )
00019 #endif
00020 //  Software Guide : BeginCommandLineArgs
00021 //     INPUTS: {FatMRISlice.png} 
00022 //     OUTPUTS: {ImageRegionIteratorOutput.png}
00023 //     20 70 210 140
00024 //  Software Guide : EndCommandLineArgs
00025 
00026 // Software Guide : BeginLatex
00027 //
00028 // \index{Iterators!speed}
00029 // The \doxygen{ImageRegionIterator} is optimized for
00030 // iteration speed and is the first choice for iterative, pixel-wise operations
00031 // when location in the image is not
00032 // important. ImageRegionIterator is the least specialized of the ITK
00033 // image iterator classes.  It implements all of the methods described in the
00034 // preceding section.
00035 //
00036 // The following example illustrates the use of
00037 // \doxygen{ImageRegionConstIterator} and ImageRegionIterator.
00038 // Most of the code constructs introduced apply to other ITK iterators as
00039 // well. This simple application crops a subregion from an image by copying
00040 // its pixel values into to a second, smaller image.
00041 //
00042 // \index{Iterators!and image regions}
00043 // \index{itk::ImageRegionIterator!example of using|(}
00044 // We begin by including the appropriate header files.
00045 // Software Guide : EndLatex
00046 
00047 #include "itkImage.h"
00048 // Software Guide : BeginCodeSnippet
00049 #include "itkImageRegionConstIterator.h"
00050 #include "itkImageRegionIterator.h"
00051 // Software Guide : EndCodeSnippet
00052 #include "itkImageFileReader.h"
00053 #include "itkImageFileWriter.h"
00054 
00055 int main( int argc, char *argv[] )
00056 {
00057   // Verify the number of parameters on the command line.
00058   if ( argc < 7 )
00059     {
00060       std::cerr << "Missing parameters. " << std::endl;
00061       std::cerr << "Usage: " << std::endl;
00062       std::cerr << argv[0]
00063                 << " inputImageFile outputImageFile startX startY sizeX sizeY"
00064                 << std::endl;
00065       return -1;
00066     }
00067 
00068 // Software Guide : BeginLatex
00069 //
00070 // Next we define a pixel type and corresponding image type. ITK iterator
00071 // classes expect the image type as their template parameter.
00072 //
00073 // Software Guide : EndLatex
00074 
00075   // Software Guide : BeginCodeSnippet
00076   const unsigned int Dimension = 2;
00077   
00078   typedef unsigned char PixelType;
00079   typedef itk::Image< PixelType, Dimension >  ImageType;
00080   
00081   typedef itk::ImageRegionConstIterator< ImageType > ConstIteratorType;
00082   typedef itk::ImageRegionIterator< ImageType>       IteratorType;
00083   // Software Guide : EndCodeSnippet
00084   
00085   typedef itk::ImageFileReader< ImageType > ReaderType;
00086   typedef itk::ImageFileWriter< ImageType > WriterType;
00087 
00088 // Software Guide : BeginLatex
00089 // 
00090 // Information about the subregion to copy is read from the command line. The
00091 // subregion is defined by an \doxygen{ImageRegion} object, with a starting
00092 // grid index and a size (Section~\ref{sec:ImageSection}).
00093 //
00094 // Software Guide : EndLatex
00095 
00096 // Software Guide : BeginCodeSnippet
00097   ImageType::RegionType inputRegion;
00098 
00099   ImageType::RegionType::IndexType inputStart;
00100   ImageType::RegionType::SizeType  size;
00101 
00102   inputStart[0] = ::atoi( argv[3] );
00103   inputStart[1] = ::atoi( argv[4] );
00104 
00105   size[0]  = ::atoi( argv[5] );
00106   size[1]  = ::atoi( argv[6] );
00107 
00108   inputRegion.SetSize( size );
00109   inputRegion.SetIndex( inputStart );
00110 // Software Guide : EndCodeSnippet
00111 
00112 
00113 // Software Guide : BeginLatex
00114 // 
00115 // The destination region in the output image is defined using the input region
00116 // size, but a different start index.  The starting index for the destination
00117 // region is the corner of the newly generated image.
00118 //
00119 // Software Guide : EndLatex
00120 
00121 // Software Guide : BeginCodeSnippet
00122   ImageType::RegionType outputRegion;
00123 
00124   ImageType::RegionType::IndexType outputStart;
00125 
00126   outputStart[0] = 0;
00127   outputStart[1] = 0;
00128 
00129   outputRegion.SetSize( size );
00130   outputRegion.SetIndex( outputStart );
00131 // Software Guide : EndCodeSnippet
00132 
00133 
00134   ReaderType::Pointer reader = ReaderType::New();
00135   reader->SetFileName( argv[1] );
00136   try
00137     {
00138     reader->Update();
00139     }
00140   catch ( itk::ExceptionObject &err)
00141     {
00142     std::cerr << "ExceptionObject caught !" << std::endl; 
00143     std::cerr << err << std::endl; 
00144     return -1;
00145     }
00146 
00147   // Check that the region is contained within the input image.
00148   if ( ! reader->GetOutput()->GetRequestedRegion().IsInside( inputRegion ) )
00149     {
00150     std::cerr << "Error" << std::endl;
00151     std::cerr << "The region " << inputRegion << "is not contained within the input image region "
00152               << reader->GetOutput()->GetRequestedRegion() << std::endl;
00153     return -1;
00154     }
00155 
00156 // Software Guide : BeginLatex
00157 //
00158 // After reading the input image and checking that the desired subregion is,
00159 // in fact, contained in the input, we allocate an output image.  It is
00160 // fundamental to set valid values to some of the basic image information 
00161 // during the copying process.  
00162 // In particular, the starting index of the output region
00163 // is now filled up with zero values and the coordinates of the physical
00164 // origin are computed as a shift from the origin of the input image. This is
00165 // quite important since it will allow us to later 
00166 // register the extracted region against the original image.
00167 //
00168 // Software Guide : EndLatex
00169 
00170 // Software Guide : BeginCodeSnippet
00171   ImageType::Pointer outputImage = ImageType::New();
00172   outputImage->SetRegions( outputRegion );
00173   const ImageType::SpacingType& spacing = reader->GetOutput()->GetSpacing();
00174   const ImageType::PointType& inputOrigin = reader->GetOutput()->GetOrigin();
00175   double   outputOrigin[ Dimension ];
00176 
00177   for(unsigned int i=0; i< Dimension; i++)
00178     {
00179     outputOrigin[i] = inputOrigin[i] + spacing[i] * inputStart[i];
00180     }
00181 
00182   outputImage->SetSpacing( spacing );
00183   outputImage->SetOrigin(  outputOrigin );
00184   outputImage->Allocate();
00185 // Software Guide : EndCodeSnippet
00186 
00187 
00188 // Software Guide : BeginLatex
00189 //
00190 // \index{Iterators!construction of} \index{Iterators!and image regions} 
00191 // The necessary images and region definitions are now in place.  All that is
00192 // left to do is to create the iterators and perform the copy.  Note that image
00193 // iterators are not accessed via smart pointers so they are light-weight
00194 // objects that are instantiated on the stack.  Also notice how the input and
00195 // output iterators are defined over the \emph{same corresponding region}.  Though the
00196 // images are different sizes, they both contain the same target subregion.
00197 //
00198 // Software Guide : EndLatex
00199 
00200 // Software Guide : BeginCodeSnippet
00201   ConstIteratorType inputIt(   reader->GetOutput(), inputRegion  );
00202   IteratorType      outputIt(  outputImage,         outputRegion );
00203 
00204   for ( inputIt.GoToBegin(), outputIt.GoToBegin(); !inputIt.IsAtEnd();
00205         ++inputIt, ++outputIt)
00206     {
00207     outputIt.Set(  inputIt.Get()  );
00208     }
00209 // Software Guide : EndCodeSnippet
00210 
00211 
00212 // Software Guide : BeginLatex
00213 //
00214 // \index{Iterators!image dimensionality}
00215 //  The \code{for} loop above is a common
00216 // construct in ITK.  The beauty of these four lines of code is that they are
00217 // equally valid for one, two, three, or even ten dimensional data, and no
00218 // knowledge of the size of the image is necessary.  Consider the ugly
00219 // alternative of ten nested \code{for} loops for traversing an image.
00220 //
00221 // Software Guide : EndLatex
00222   
00223   WriterType::Pointer writer = WriterType::New();
00224   writer->SetFileName( argv[2] );
00225   writer->SetInput( outputImage );
00226 
00227   try
00228     {
00229     writer->Update();
00230     }
00231   catch ( itk::ExceptionObject &err)
00232     {
00233     std::cerr << "ExceptionObject caught !" << std::endl; 
00234     std::cerr << err << std::endl; 
00235     return -1;   
00236     }
00237 
00238 // Software Guide : BeginLatex
00239 //
00240 // Let's run this example on the image \code{FatMRISlice.png} found
00241 // in \code{Examples/Data}.  The command line arguments specify the
00242 // input and output file names, then the $x$, $y$ origin and the $x$, $y$ size
00243 // of the cropped subregion.
00244 //
00245 // \small
00246 // \begin{verbatim}
00247 // ImageRegionIterator FatMRISlice.png ImageRegionIteratorOutput.png 20 70 210 140
00248 // \end{verbatim}
00249 // \normalsize
00250 //
00251 // The output is the cropped subregion shown in
00252 // Figure~\ref{fig:ImageRegionIteratorOutput}.
00253 //
00254 // \begin{figure}
00255 // \centering
00256 // \includegraphics[width=0.4\textwidth]{FatMRISlice.eps}
00257 // \includegraphics[width=0.3\textwidth]{ImageRegionIteratorOutput.eps}
00258 // \itkcaption[Copying an image subregion using ImageRegionIterator]{Cropping a
00259 // region from an image.  The original image is shown at left.  The image on
00260 // the right is the result of applying the ImageRegionIterator example code.}
00261 // \protect\label{fig:ImageRegionIteratorOutput}
00262 // \end{figure}
00263 // 
00264 // \index{itk::ImageRegionIterator!example of using|)}
00265 //
00266 // Software Guide : EndLatex
00267 
00268   return 0;
00269 }
00270 
00271 

Generated at Sun Sep 23 11:49:07 2007 for ITK by doxygen 1.5.1 written by Dimitri van Heesch, © 1997-2000