[Insight-users] RGBAPixel Images

Matthew Gibb matthew.gibb at dtc.ox.ac.uk
Tue Jan 18 18:50:17 EST 2011


Hello Insight mailing list!

I sent a copy of this email a couple of days ago, but I think it's been
rejected because I attached a rather large image to it, so apologies if two
of these come through :-)

I think have either found a bug in ITK, or I am doing something silly.
Either way, I'm really confused and I'd really appreciate some advice...

I need to downsample two sets images of heart data and register them
together, but one set consists of RGB bmps, and the other of RGBA bmps. A
example of one of the original RGBA images can be found here:

http://dl.dropbox.com/u/706648/0590.bmp

and an example of the output of an ITK program that has the problems I
discuss can be found here:

http://dl.dropbox.com/u/706648/0590_ITK_output.bmp

When I try and read RGBA bmps with a reader of pixel type RGB, the output
looks very different, because (I guess) the actual alpha channel has been
interpreted as red, the actual red as green, and the actual green as blue
(the original blue channel seems to be discarded).

Furthermore, when I read the RGBA bmp with a reader of pixel type RGBA, then
explicitly populate new images of type itkImage< unsigned char, 2 > using
either the number getter GetNthComponent() or the named getters GetRed(),
GetGreen(), GetBlue() and GetAlpha(), as in the code below, the same problem
is manifest:

// ConvertRGBAtoRGB.cxx
#include "itkRGBPixel.h"
#include "itkRGBAPixel.h"
#include "IOHelpers.hpp"


using namespace std;

int main( int argc, char * argv[] )
{
  if( argc != 3 ) {
    cerr << "Usage: " << endl;
    cerr << argv[0] << "  inputImageFile outputDir" << endl;
    exit(EXIT_FAILURE);
  }

  string inputImageFile(argv[1]), outputDir(argv[2]);

  const     unsigned int   Dimension = 2;

  typedef itk::RGBAPixel< unsigned char > InputPixelType;
  typedef itk::RGBPixel< unsigned char >  OutputPixelType;

  typedef itk::Image< InputPixelType,  Dimension > InputImageType;
  typedef itk::Image< OutputPixelType, Dimension > OutputImageType;

  // Read RGBA image from file
  InputImageType::Pointer inputImage = readImage< InputImageType >(
inputImageFile );

  // Build RGB image from RGBA image
  // Initialise output image
  OutputImageType::RegionType region;
  region.SetSize( inputImage->GetLargestPossibleRegion().GetSize() );

  OutputImageType::Pointer outputImage = OutputImageType::New();
  outputImage->SetRegions( region );
  outputImage->CopyInformation( inputImage );
  outputImage->Allocate();

  // copy values from input image
  itk::ImageRegionConstIterator< InputImageType > cit(inputImage, region);
  itk::ImageRegionIterator< OutputImageType >     it(outputImage, region);
  for (cit.GoToBegin(), it.GoToBegin(); !it.IsAtEnd(); ++cit, ++it ) {
    it.Value().SetRed( cit.Value().GetRed() );
    it.Value().SetBlue( cit.Value().GetBlue() );
    it.Value().SetGreen( cit.Value().GetGreen() );

  }

  // write RGB image to file
  writeImage< OutputImageType >( outputImage, outputDir + "RGB.bmp" );

  // Write individual channels by name to file.
  #define writeChannelByName(Colour) \
  { \
   /* create image from named channel of input image */ \
    typedef itk::Image< unsigned char > ChannelType; \
    ChannelType::Pointer channel = ChannelType::New(); \
    channel->SetRegions( region ); \
    channel->CopyInformation( inputImage ); \
    channel->Allocate(); \
    \
    itk::ImageRegionConstIterator< InputImageType > cit(inputImage, region);
\
    itk::ImageRegionIterator< ChannelType > it(channel, region); \
    \
    for (cit.GoToBegin(), it.GoToBegin(); !it.IsAtEnd(); ++cit, ++it ) { \
      it.Set( cit.Value().Get ## Colour() ); \
    } \
    \
   /* write channel */ \
    writeImage< ChannelType >(channel, outputDir + #Colour + ".bmp" ); \
  }

  writeChannelByName(Red);
  writeChannelByName(Green);
  writeChannelByName(Blue);
  writeChannelByName(Alpha);

  return EXIT_SUCCESS;
}

// IOHelpers.hpp
#ifndef IO_HELPERS_HPP_
#define IO_HELPERS_HPP_

#include "itkImage.h"
#include "itkImageFileReader.h"
#include "itkImageFileWriter.h"
#include "itkTransformFileWriter.h"

using namespace std;

// Reader helper
template<typename ImageType>
typename ImageType::Pointer readImage(const string& fileName) {
  typedef itk::ImageFileReader< ImageType > ReaderType;
  typename ReaderType::Pointer reader = ReaderType::New();
  reader->SetFileName( fileName.c_str() );
   try {
   reader->Update();
  }
catch( itk::ExceptionObject & err ) {
    cerr << "ExceptionObject caught !" << endl;
    cerr << err << endl;
exit(EXIT_FAILURE);
}
   return reader->GetOutput();
}

// Writer helpers
// const Data
template<typename WriterType, typename DataType>
void writeData(const typename DataType::ConstPointer data, const string&
fileName) {
  typename WriterType::Pointer writer = WriterType::New();
 writer->SetInput( data );
  // writer->AddTransform( anotherTransform ); // only applies to writing
transforms

  writer->SetFileName( fileName.c_str() );
   try {
   writer->Update();
  }
catch( itk::ExceptionObject & err ) {
    cerr << "ExceptionObject caught !" << endl;
    cerr << err << endl;
exit(EXIT_FAILURE);
}
}

// Data
template<typename WriterType, typename DataType>
void writeData(const typename DataType::Pointer data, const string&
fileName) {
  writeData< WriterType, DataType >( (typename DataType::ConstPointer) data,
fileName);
}

// Const Image
template<typename ImageType>
void writeImage(const typename ImageType::ConstPointer image, const string&
fileName) {
  writeData< itk::ImageFileWriter< ImageType >, ImageType >( image, fileName
);
}

// Image
template<typename ImageType>
void writeImage(const typename ImageType::Pointer image, const string&
fileName) {
  writeData< itk::ImageFileWriter< ImageType >, ImageType >( image, fileName
);
}

#endif

It seems like the internal RGBA image pointer is one char behind where it
should be or something...?

I am using the current 3.20 release branch from the Git repo on a Penryn
MacBook Pro running Snow Leopard, but I get the same problem on other linux
boxes.

Unrelatedly, I'm always keen for feedback on coding style, conformance,
anything really, so please don't hold back!

Many thanks,

Matt Gibb.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://www.itk.org/pipermail/insight-users/attachments/20110118/987e5704/attachment.htm>


More information about the Insight-users mailing list