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

Iterators/ImageRegionIteratorWithIndex.cxx

For a complete description of the ITK Image Iterators and their API, please see the Iterators chapter in the ITK Software Guide. The ITK Software Guide is available in print and as a free .pdf download from http://www.itk.org.

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: ImageRegionIteratorWithIndex.cxx,v $
00005   Language:  C++
00006   Date:      $Date: 2005/02/08 03:59:00 $
00007   Version:   $Revision: 1.25 $
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: {VisibleWomanEyeSlice.png} 
00022 //     OUTPUTS: {ImageRegionIteratorWithIndexOutput.png}
00023 //  Software Guide : EndCommandLineArgs
00024 
00025 // Software Guide : BeginLatex
00026 //
00027 // \index{Iterators!speed}
00028 // The ``WithIndex'' family of iterators was designed for algorithms that
00029 // use both the value and the location of image pixels in calculations.  Unlike
00030 // \doxygen{ImageRegionIterator}, which calculates an index only when
00031 // asked for, \doxygen{ImageRegionIteratorWithIndex} maintains its
00032 // index location as a member variable that is updated during the increment or
00033 // decrement process. Iteration speed is penalized, but the index queries are
00034 // more efficient.
00035 //
00036 // \index{itk::ImageRegionIteratorWithIndex!example of using|(}
00037 //
00038 // The following example illustrates the use of
00039 // ImageRegionIteratorWithIndex.  The algorithm mirrors
00040 // a 2D image across its $x$-axis (see \doxygen{FlipImageFilter} for an ND
00041 // version).  The algorithm makes extensive use of the \code{GetIndex()}
00042 // method.
00043 //
00044 // We start by including the proper header file.
00045 //
00046 // Software Guide : EndLatex
00047 
00048 #include "itkImage.h"
00049 #include "itkRGBPixel.h"
00050 // Software Guide : BeginCodeSnippet
00051 #include "itkImageRegionIteratorWithIndex.h"
00052 // Software Guide : EndCodeSnippet
00053 #include "itkImageFileReader.h"
00054 #include "itkImageFileWriter.h"
00055 
00056 int main( int argc, char *argv[] )
00057 {
00058   // Verify the number of parameters on the command line.
00059   if ( argc < 3 )
00060     {
00061     std::cerr << "Missing parameters. " << std::endl;
00062     std::cerr << "Usage: " << std::endl;
00063     std::cerr << argv[0]
00064               << " inputImageFile outputImageFile"
00065               << std::endl;
00066     return -1;
00067     }
00068 
00069 // Software Guide : BeginLatex
00070 //
00071 // For this example, we will use an RGB pixel type so that we can process color
00072 // images. Like most other ITK image iterator,
00073 // ImageRegionIteratorWithIndex class expects the image type as its
00074 // single template parameter.
00075 //
00076 // Software Guide : EndLatex
00077 
00078 // Software Guide : BeginCodeSnippet
00079   const unsigned int Dimension = 2;
00080   
00081   typedef itk::RGBPixel< unsigned char > RGBPixelType;
00082   typedef itk::Image< RGBPixelType, Dimension >  ImageType;
00083   
00084   typedef itk::ImageRegionIteratorWithIndex< ImageType > IteratorType;
00085 // Software Guide : EndCodeSnippet
00086   
00087   typedef itk::ImageFileReader< ImageType > ReaderType;
00088   typedef itk::ImageFileWriter< ImageType > WriterType;
00089 
00090   ImageType::ConstPointer inputImage;
00091   ReaderType::Pointer reader = ReaderType::New();
00092   reader->SetFileName( argv[1] );
00093   try
00094     {
00095     reader->Update();
00096     inputImage = reader->GetOutput();
00097     }
00098   catch ( itk::ExceptionObject &err)
00099     {
00100     std::cout << "ExceptionObject caught !" << std::endl; 
00101     std::cout << err << std::endl; 
00102     return -1;
00103     }
00104 
00105 // Software Guide : BeginLatex
00106 //
00107 // An \code{ImageType} smart pointer called \code{inputImage} points to the
00108 // output of the image reader.  After updating the image reader, we can
00109 // allocate an output image of the same size, spacing, and origin as the
00110 // input image.
00111 //
00112 // Software Guide : EndLatex
00113 
00114 // Software Guide : BeginCodeSnippet
00115   ImageType::Pointer outputImage = ImageType::New();
00116   outputImage->SetRegions( inputImage->GetRequestedRegion() );
00117   outputImage->CopyInformation( inputImage );
00118   outputImage->Allocate();
00119 // Software Guide : EndCodeSnippet
00120 
00121 // Software Guide : BeginLatex
00122 //
00123 // Next we create the iterator that walks the output image.  This algorithm
00124 // requires no iterator for the input image.
00125 //
00126 // Software Guide : EndLatex
00127 
00128 // Software Guide : BeginCodeSnippet
00129   IteratorType outputIt( outputImage, outputImage->GetRequestedRegion() );
00130 // Software Guide : EndCodeSnippet
00131 
00132 // Software Guide: BeginLatex
00133 //
00134 // This axis flipping algorithm works by iterating through the output image,
00135 // querying the iterator for its index, and copying the value from the input
00136 // at an index mirrored across the $x$-axis.
00137 //
00138 // Software Guide : EndLatex
00139 
00140 // Software Guide : BeginCodeSnippet
00141   ImageType::IndexType requestedIndex =
00142                 outputImage->GetRequestedRegion().GetIndex();
00143   ImageType::SizeType requestedSize =
00144                 outputImage->GetRequestedRegion().GetSize();
00145 
00146   for ( outputIt.GoToBegin(); !outputIt.IsAtEnd(); ++outputIt)
00147     {
00148     ImageType::IndexType idx = outputIt.GetIndex();
00149     idx[0] =  requestedIndex[0] + requestedSize[0] - 1 - idx[0];
00150     outputIt.Set( inputImage->GetPixel(idx) );
00151     }
00152 // Software Guide : EndCodeSnippet
00153 
00154   WriterType::Pointer writer = WriterType::New();
00155   writer->SetFileName( argv[2] );
00156   writer->SetInput(outputImage);
00157   try
00158     {
00159     writer->Update();
00160     }
00161   catch ( itk::ExceptionObject &err)
00162     {
00163     std::cout << "ExceptionObject caught !" << std::endl; 
00164     std::cout << err << std::endl; 
00165     return -1;   
00166 }
00167 
00168 // Software Guide : BeginLatex
00169 //
00170 // Let's run this example on the image \code{VisibleWomanEyeSlice.png} found in
00171 // the \code{Examples/Data} directory.
00172 // Figure~\ref{fig:ImageRegionIteratorWithIndexExample} shows how the original
00173 // image has been mirrored across its $x$-axis in the output.
00174 //
00175 // \begin{figure} \center
00176 // \includegraphics[width=0.44\textwidth]{VisibleWomanEyeSlice.eps}
00177 // \includegraphics[width=0.44\textwidth]{ImageRegionIteratorWithIndexOutput.eps}
00178 // \itkcaption[Using the ImageRegionIteratorWithIndex]{Results of using
00179 // ImageRegionIteratorWithIndex to mirror an image across an axis. The original
00180 // image is shown at left.  The mirrored output is shown at right.}
00181 // \label{fig:ImageRegionIteratorWithIndexExample}
00182 // \end{figure}
00183 //
00184 // \index{itk::ImageRegionIteratorWithIndex!example of using|)}
00185 //
00186 // Software Guide : EndLatex
00187   
00188   return 0;
00189 }

Generated at Tue Jul 29 19:03:31 2008 for ITK by doxygen 1.5.1 written by Dimitri van Heesch, © 1997-2000