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

DataRepresentation/Image/RGBImage.cxx

00001 /*=========================================================================
00002 
00003   Program:   Insight Segmentation & Registration Toolkit
00004   Module:    $RCSfile: RGBImage.cxx,v $
00005   Language:  C++
00006   Date:      $Date: 2005/02/08 03:51:53 $
00007   Version:   $Revision: 1.19 $
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 #include "itkImage.h"
00022 #include "itkImageFileReader.h"
00023 
00024 //  Software Guide : BeginLatex
00025 //
00026 //  Thanks to the flexibility offered by the
00027 //  \href{http://www.boost.org/more/generic_programming.html}{Generic
00028 //  Programming} style on which ITK is based, it is possible to
00029 //  instantiate images of arbitrary pixel type.  The following example
00030 //  illustrates how a color image with RGB pixels can be defined.
00031 //
00032 //  A class intended to support the RGB pixel type is available in ITK.  You
00033 //  could also define your own pixel class and use it to instantiate a
00034 //  custom image type. In order to use the \doxygen{RGBPixel} class, it is
00035 //  necessary to include its header file.
00036 //
00037 //  \index{itk::RGBPixel}
00038 //  \index{itk::RGBPixel!Image}
00039 //  \index{itk::RGBPixel!header}
00040 //
00041 //  Software Guide : EndLatex 
00042 
00043 // Software Guide : BeginCodeSnippet
00044 #include "itkRGBPixel.h"
00045 // Software Guide : EndCodeSnippet
00046 
00047 
00048 int main( int , char * argv[] )
00049 {
00050   // Software Guide : BeginLatex
00051   //
00052   // The RGB pixel class is templated over a type used to represent each one
00053   // of the red, green and blue pixel components. A typical instantiation of the
00054   // templated class is as follows.
00055   //
00056   //  \index{itk::RGBPixel!Instantiation}
00057   //
00058   // Software Guide : EndLatex 
00059 
00060   // Software Guide : BeginCodeSnippet
00061   typedef itk::RGBPixel< unsigned char >    PixelType;
00062   // Software Guide : EndCodeSnippet
00063 
00064 
00065   // Software Guide : BeginLatex
00066   //
00067   // The type is then used as the pixel template parameter of the image.
00068   //
00069   // Software Guide : EndLatex 
00070 
00071   // Software Guide : BeginCodeSnippet
00072   typedef itk::Image< PixelType, 3 >   ImageType;
00073   // Software Guide : EndCodeSnippet
00074 
00075 
00076   // Software Guide : BeginLatex
00077   //
00078   // The image type can be used to instantiate other filter, for example,
00079   // an \doxygen{ImageFileReader} object that will read the image from a
00080   // file.
00081   //
00082   // \index{itk::ImageFileReader!RGB Image}
00083   //
00084   // Software Guide : EndLatex 
00085 
00086   // Software Guide : BeginCodeSnippet
00087   typedef itk::ImageFileReader< ImageType >  ReaderType;
00088   // Software Guide : EndCodeSnippet
00089 
00090 
00091   ReaderType::Pointer reader = ReaderType::New();
00092   const char * filename = argv[1];
00093   reader->SetFileName( filename );
00094   reader->Update();
00095 
00096   ImageType::Pointer image = reader->GetOutput();
00097 
00098   ImageType::IndexType pixelIndex;
00099 
00100   pixelIndex[0] = 25;  
00101   pixelIndex[1] = 35;  
00102   pixelIndex[2] =  0;  
00103 
00104 
00105   // Software Guide : BeginLatex
00106   //
00107   // Access to the color components of the pixels can now be performed using
00108   // the methods provided by the RGBPixel class.
00109   //
00110   // \index{itk::Image!GetPixel()}
00111   // \index{itk::RGBPixel!GetRed()}
00112   // \index{itk::RGBPixel!GetGreen()}
00113   // \index{itk::RGBPixel!GetBlue()}
00114   //
00115   // Software Guide : EndLatex 
00116 
00117   // Software Guide : BeginCodeSnippet
00118   PixelType onePixel = image->GetPixel( pixelIndex );
00119   
00120   PixelType::ValueType red   = onePixel.GetRed();
00121   PixelType::ValueType green = onePixel.GetGreen();
00122   PixelType::ValueType blue  = onePixel.GetBlue();
00123   // Software Guide : EndCodeSnippet
00124 
00125   std::cout << "Pixel values from GetRed,GetGreen,GetBlue:" << std::endl;
00126   std::cout << "Red = "
00127             << itk::NumericTraits<PixelType::ValueType>::PrintType(red)
00128             << std::endl;
00129   std::cout << "Green = "
00130             << itk::NumericTraits<PixelType::ValueType>::PrintType(green)
00131             << std::endl;
00132   std::cout << "Blue = "
00133             << itk::NumericTraits<PixelType::ValueType>::PrintType(blue)
00134             << std::endl;
00135 
00136   // Software Guide : BeginLatex
00137   //
00138   // The subindex notation can also be used since the \doxygen{RGBPixel} inherits the
00139   // \code{[]} operator from the \doxygen{FixedArray} class.
00140   //
00141   // Software Guide : EndLatex 
00142 
00143   // Software Guide : BeginCodeSnippet
00144   red   = onePixel[0];  // extract Red   component
00145   green = onePixel[1];  // extract Green component
00146   blue  = onePixel[2];  // extract Blue  component
00147 
00148   std::cout << "Pixel values:" << std::endl;
00149   std::cout << "Red = "
00150             << itk::NumericTraits<PixelType::ValueType>::PrintType(red)
00151             << std::endl;
00152   std::cout << "Green = "
00153             << itk::NumericTraits<PixelType::ValueType>::PrintType(green)
00154             << std::endl;
00155   std::cout << "Blue = "
00156             << itk::NumericTraits<PixelType::ValueType>::PrintType(blue)
00157             << std::endl;
00158 
00159   // Software Guide : EndCodeSnippet
00160 
00161   // Lets repeat that both \code{SetPixel()} and \code{GetPixel()} are
00162   // inefficient and should only be used for debugging purposes or for
00163   // implementing interactions with a graphical user interface such as
00164   // querying pixel value by clicking with the mouse.
00165   //
00166  
00167   return 0;
00168 }
00169 

Generated at Wed Nov 5 20:13:04 2008 for ITK by doxygen 1.5.1 written by Dimitri van Heesch, © 1997-2000