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

Examples/Iterators/ImageLinearIteratorWithIndex.cxx

MORE INFORMATION
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: ImageLinearIteratorWithIndex.cxx,v $
00005   Language:  C++
00006   Date:      $Date: 2008-01-23 16:43:05 $
00007   Version:   $Revision: 1.23 $
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 
00021 // Software Guide : BeginLatex
00022 //
00023 // The \doxygen{ImageLinearIteratorWithIndex} is designed for line-by-line
00024 // processing of an image.  It walks a linear path along a selected image
00025 // direction parallel to one of the coordinate axes of the image. This
00026 // iterator conceptually breaks an image into a set of parallel lines
00027 // that span the selected image dimension.
00028 //
00029 // \index{Iterators!and image lines}
00030 //
00031 // Like all image iterators, movement of the
00032 // ImageLinearIteratorWithIndex is constrained within an
00033 // image region $R$.  The line $\ell$ through which the iterator moves is
00034 // defined by selecting a direction and an origin.   The line $\ell$
00035 // extends from the origin to the upper boundary of $R$. The origin can be
00036 // moved to any position along the lower boundary of $R$. 
00037 //    
00038 // Several additional methods are defined for this iterator to control movement
00039 // of the iterator along the line $\ell$ and movement of the origin of $\ell$.
00040 //
00041 // %Might need a figure here to describe this iterator.
00042 //
00043 // \begin{itemize}
00044 //
00045 // \index{itk::ImageLinearIteratorWithIndex!NextLine()}
00046 //
00047 // \item \textbf{\code{NextLine()}} Moves the iterator to the beginning pixel
00048 // location of the next line in the image.  The origin of the next line is
00049 // determined by incrementing the current origin along the fastest increasing
00050 // dimension of the subspace of the image that excludes the selected dimension.
00051 //
00052 //
00053 // \index{itk::ImageLinearIteratorWithIndex!PreviousLine()}
00054 //
00055 // \item \textbf{\code{PreviousLine()}} Moves the iterator to the \emph{last valid
00056 // pixel location} in the previous line. The origin of the previous line is
00057 // determined by decrementing the current origin along the fastest increasing
00058 // dimension of the subspace of the image that excludes the selected dimension.
00059 //
00060 // \index{itk::ImageLinearIteratorWithIndex!GoToBeginOfLine()}
00061 //
00062 // \item \textbf{\code{GoToBeginOfLine()}} Moves the iterator to the beginning
00063 // pixel of the current line.
00064 //
00065 // \index{itk::ImageLinearIteratorWithIndex!GoToEndOfLine()}
00066 //
00067 // \item \textbf{\code{GoToEndOfLine()}}  Move the iterator to 
00068 // \emph{one past} the last valid pixel of the current line.
00069 //
00070 // \index{itk::ImageLinearIteratorWithIndex!GoToReverseBeginOfLine()}
00071 //
00072 // \item \textbf{\code{GoToReverseBeginOfLine()}}  Move the iterator
00073 // to \emph{the last valid pixel} of the current line.
00074 //
00075 // \index{itk::ImageLinearIteratorWithIndex!IsAtReverseEndOfLine()}
00076 //
00077 // \item \textbf{\code{IsAtReverseEndOfLine()}} 
00078 // Returns true if the iterator points
00079 // to \emph{one position before} the beginning pixel of the current line.
00080 //
00081 // \index{itk::ImageLinearIteratorWithIndex!IsAtEndOfLine()}
00082 //
00083 // \item \textbf{\code{IsAtEndOfLine()}}  
00084 // Returns true if the iterator points to
00085 // \emph{one position past} the last valid pixel of the current line.
00086 // \end{itemize}
00087 //
00088 // The following code example shows how to use the
00089 // ImageLinearIteratorWithIndex.  It implements the same algorithm as
00090 // in the previous example, flipping an image across its $x$-axis.  Two line
00091 // iterators are iterated in opposite directions across the $x$-axis.  
00092 // After each line is traversed, the iterator origins are stepped along 
00093 // the $y$-axis to the
00094 // next line.
00095 //
00096 // \index{itk::ImageLinearIteratorWithIndex!example of using|(}
00097 //
00098 // Headers for both the const and non-const versions are needed.
00099 //
00100 // Software Guide : EndLatex
00101 
00102 #include "itkImage.h"
00103 #include "itkRGBPixel.h"
00104 // Software Guide : BeginCodeSnippet
00105 #include "itkImageLinearConstIteratorWithIndex.h"
00106 #include "itkImageLinearIteratorWithIndex.h"
00107 // Software Guide : EndCodeSnippet
00108 #include "itkImageFileReader.h"
00109 #include "itkImageFileWriter.h"
00110 
00111 int main( int argc, char *argv[] )
00112 {
00113   // Verify the number of parameters on the command line.
00114   if ( argc < 3 )
00115     {
00116     std::cerr << "Missing parameters. " << std::endl;
00117     std::cerr << "Usage: " << std::endl;
00118     std::cerr << argv[0]
00119               << " inputImageFile outputImageFile"
00120               << std::endl;
00121     return -1;
00122     }
00123 
00124 // Software Guide : BeginLatex
00125 //
00126 // The RGB image and pixel types are defined as in the previous example.  The 
00127 // ImageLinearIteratorWithIndex class and its const version each have
00128 // single template parameters, the image type.
00129 //
00130 // Software Guide : EndLatex
00131 
00132   const unsigned int Dimension = 2;
00133   
00134   typedef itk::RGBPixel< unsigned char > RGBPixelType;
00135   typedef itk::Image< RGBPixelType, Dimension >  ImageType;
00136 
00137 // Software Guide : BeginCodeSnippet
00138   typedef itk::ImageLinearIteratorWithIndex< ImageType >       IteratorType;
00139   typedef itk::ImageLinearConstIteratorWithIndex< ImageType >  ConstIteratorType;
00140 // Software Guide : EndCodeSnippet
00141   
00142   typedef itk::ImageFileReader< ImageType > ReaderType;
00143   typedef itk::ImageFileWriter< ImageType > WriterType;
00144 
00145   ImageType::ConstPointer inputImage;
00146   ReaderType::Pointer reader = ReaderType::New();
00147   reader->SetFileName( argv[1] );
00148   try
00149     {
00150     reader->Update();
00151     inputImage = reader->GetOutput();
00152     }
00153   catch ( itk::ExceptionObject &err)
00154     {
00155     std::cout << "ExceptionObject caught a !" << std::endl; 
00156     std::cout << err << std::endl; 
00157     return -1;
00158     }
00159 
00160 // Software Guide : BeginLatex
00161 //
00162 // After reading the input image, we allocate an output image that of the same
00163 // size, spacing, and origin.
00164 //
00165 // Software Guide : EndLatex
00166 
00167 // Software Guide : BeginCodeSnippet
00168   ImageType::Pointer outputImage = ImageType::New();
00169   outputImage->SetRegions( inputImage->GetRequestedRegion() );
00170   outputImage->CopyInformation( inputImage );
00171   outputImage->Allocate();
00172 // Software Guide : EndCodeSnippet
00173 
00174 // Software Guide : BeginLatex
00175 //
00176 // Next we create the two iterators.  The const iterator walks the input image,
00177 // and the non-const iterator walks the output image.  The iterators are
00178 // initialized over the same region.  The direction of iteration is set to 0,
00179 // the $x$ dimension.
00180 //
00181 // Software Guide : EndLatex
00182 
00183 // Software Guide : BeginCodeSnippet
00184   ConstIteratorType inputIt( inputImage, inputImage->GetRequestedRegion() );
00185   IteratorType outputIt( outputImage, inputImage->GetRequestedRegion() );
00186 
00187   inputIt.SetDirection(0);
00188   outputIt.SetDirection(0);
00189 // Software Guide : EndCodeSnippet
00190 
00191 // Software Guide: BeginLatex
00192 //
00193 // Each line in the input is copied to the output.  The input iterator moves
00194 // forward across columns while the output iterator moves backwards.
00195 //
00196 // Software Guide : EndLatex
00197 
00198 // Software Guide : BeginCodeSnippet
00199   for ( inputIt.GoToBegin(),  outputIt.GoToBegin(); ! inputIt.IsAtEnd();
00200         outputIt.NextLine(),  inputIt.NextLine())
00201     {
00202     inputIt.GoToBeginOfLine();
00203     outputIt.GoToEndOfLine();
00204     while ( ! inputIt.IsAtEndOfLine() )
00205       {
00206       --outputIt;
00207       outputIt.Set( inputIt.Get() );
00208       ++inputIt;
00209       }
00210     }
00211 // Software Guide : EndCodeSnippet
00212 
00213   WriterType::Pointer writer = WriterType::New();
00214   writer->SetFileName( argv[2] );
00215   writer->SetInput(outputImage);
00216   try
00217     {
00218     writer->Update();
00219     }
00220   catch ( itk::ExceptionObject &err)
00221     {
00222     std::cout << "ExceptionObject caught !" << std::endl; 
00223     std::cout << err << std::endl; 
00224     return -1;   
00225     }
00226 
00227 // Software Guide : BeginLatex
00228 //
00229 // Running this example on \code{VisibleWomanEyeSlice.png} produces
00230 // the same output image shown in
00231 // Figure~\ref{fig:ImageRegionIteratorWithIndexExample}.
00232 //
00233 // \index{itk::ImageLinearIteratorWithIndex!example of using|)}
00234 // Software Guide : EndLatex
00235   
00236   return 0;
00237 }

Generated at Mon Apr 14 11:33:21 2008 for ITK by doxygen 1.5.1 written by Dimitri van Heesch, © 1997-2000