[Insight-users] how to read images with unknown number of channels at build time ?

Karthik Krishnan Karthik.Krishnan at kitware.com
Mon Mar 6 15:13:33 EST 2006


Gaetan Lehmann wrote:

>Karthik, Kent, 
>
>Thanks for your answer. I think VectorImage is what I need.
>However I'm not sure to understand the problem of pixel organization. The 
>reader with ImageVector as template argument does not work properly with all 
>file types ?
>  
>
It should work for any multi-component image. It certainly does for 
Nrrd, Analyze and MetaImage and PNG color images since there are tests 
for those, in the DTI, and Statistics folders. I haven't tested it on 
every concievable file format. If it does not work on any, please let me 
know.

---------------
Here's an exerpt from the doxygen of this class that might help augment 
documentation:

This class differs from Image 
<http://www.itk.org/Insight/Doxygen/html/classitk_1_1Image.html> in that 
it is intended to represent multiple images. Each pixel represents /k/ 
measurements, each of datatype /TPixel/. The memory organization of the 
resulting image is as follows: ... Pi0 Pi1 Pi2 Pi3 P(i+1)0 P(i+1)1 
P(i+1)2 P(i+1)3 P(i+2)0 ... where Pi0 represents the 0th measurement of 
the pixel at index i. [For instance for a 3 component color image, it 
would be   ...Rk Gk Bk R(k+1) G(k+1) B(k+1) ......]

Conceptually, a |VectorImage< double, 3 >| is the same as a |Image 
<http://www.itk.org/Insight/Doxygen/html/classitk_1_1Image.html>< 
VariableLengthVector< double >, 3 >|. The difference lies in the memory 
organization. The latter results in a fragmented organization with each 
location in the Image 
<http://www.itk.org/Insight/Doxygen/html/classitk_1_1Image.html> holding 
a pointer to an |VariableLengthVector 
<http://www.itk.org/Insight/Doxygen/html/classitk_1_1VariableLengthVector.html>| 
holding the actual pixel. The former stores the /k/ pixels instead of a 
pointer reference, which apart from avoiding fragmentation of memory 
also avoids storing a 8 bytes of pointer reference for each pixel. The 
parameter /k/ can be set using |SetVectorLength|.

The API of the class is such that it returns a pixeltype 
VariableLengthVector< double > when queried, with the data internally 
pointing to the buffer. (the container does not manage the memory). 
Similarly SetPixel calls can be made with VariableLengthVector< double >.

-----------

>Is there some doc about ImageVector ? There is nothing about it in the user 
>guide :-/
>  
>
Ok.. It needs to be added to the guide.

>Gaetan
>
>
>On Friday 03 March 2006 23:41, Karthik Krishnan wrote:
>  
>
>>Just to clarify, its the up to you to choose a representation for the
>>image (with run time variable number of channels) being read in.
>>
>>For instance given a diffusion weighted dataset such as (the header):
>>ftp://itk.org/pub/namic/DTI/Data/dwi.nhdr, or its meta-image equivalent:
>>ftp://itk.org/pub/namic/DTI/Data/xt_dwi.mhd,
>>
>>You could read such an image as a
>>1. 3D image with 31 channels -  itk::VectorImage< short, 3 >
>>2. 31 3D images   - 31 itk::Image< short, 3 >
>>3. 4D image with dims[3] = 31 - itk::Image< short, 4 >
>>
>>I assumed you wanted the first, which may not be the case.
>>
>>Remember that the organization of pixels in the file will affect what
>>you read (assuming you were reading naively at one shot..)
>>
>>1. Readers instantiated over images of kind 1 assume an order.. R G B R
>>G B R G B ...
>>2. Readers instantiated over images of kind 2 and 3 assume ..     R R R
>>...G G G... B B B ..
>>
>>Of course, you could always do your IO tricks to reorganize memory
>>locations..
>>
>>-karthik
>>
>>Kent Williams wrote:
>>    
>>
>>>It is better, on the whole, to know what you're reading and then
>>>respond accordingly.
>>>
>>>Presently there is no itk::FileSniffer, as it were.
>>>itk::ImageFileReader does try to load whatever path you ask it to and
>>>do something reasonable with the data. But there are a couple of
>>>'gotchas' if you want to keep things completely general:
>>>
>>>1. File formats where images of dimension > 3 are stored as
>>>two-dimensional slices. Then, itk::ImageSeriesReader is the tool of
>>>choice, instead of itk::ImageFileReader.
>>>2. Cases such as yours, where some attribute of the file on disk
>>>affects how you will process it.
>>>
>>>The workaround is to use an itk:ImageFileReader instance to sneak a
>>>look at the files attributes, before actually reading the image. Below
>>>is a function that illustrates how to snoop around for information
>>>about the file you're trying to load.  This feels to me like breaking
>>>encapsulation, but there currently isn't any other way.
>>>To get even more nefarious, you can call
>>>itk::ImageIOBase::GetNameOfClass() to find out which file reader
>>>claims to be able to read the file.  We use this in order to handle
>>>DICOM images specially. Our application (which didn't originally use
>>>ITK Image I/O) had the convention that when you give the name of a
>>>slice of a DICOM file, then you want to read all the slices in the
>>>containing directory that are members of the same DICOM series.
>>>
>>>int FileInfo(const std::string &path)
>>>{
>>> typedef itk::Image<char,4> TestImageType; // pixel type doesn't
>>>matter for current purpose
>>> typedef itk::ImageFileReader<TestImageType> TestFileReaderType; //
>>>reader for testing a file
>>> TestFileReaderType::Pointer onefileReader = TestFileReaderType::New();
>>>  onefileReader->SetFileName(path.c_str());
>>> try
>>>   {
>>>   onefileReader->GenerateOutputInformation();
>>>   }
>>> catch(itk::ExceptionObject &excp)
>>>   {
>>>   return -1;
>>>   }
>>>   // grab the ImageIO instance for the reader
>>> itk::ImageIOBase *imageIO = onefileReader->GetImageIO();
>>>   unsigned int NumberOfDimensions =  imageIO->GetNumberOfDimensions()
>>>std::endl;
>>>   unsigned dims[32];   // almost always no more than 4 dims, but ...
>>>   unsigned origin[32];
>>>   double spacing[32];
>>>   std::vector<double> directions[32];
>>>   for(unsigned i = 0; i < NumberOfDimensions && i < 32; i++)
>>>     {
>>>     dims[i] = imageIO->GetDimensions(i);
>>>     origin[i] = imageIO->GetOrigin(i);
>>>     spacing[i] = imageIO->GetSpacing(i);
>>>     directions[i] = imageIO->GetDirection(i);
>>>     }
>>>  // PixelType is SCALAR, RGB, RGBA, VECTOR, COVARIANTVECTOR, POINT,
>>>INDEX
>>>  itk::ImageIOBase::PixelType pixelType = imageIO->GetPixelType();
>>>  // IOComponentType is UCHAR, CHAR, USHORT, SHORT, UINT, INT, ULONG,
>>>LONG, FLOAT, DOUBLE
>>>  itk::ImageIOBase::IOComponentType componentType =
>>>imageIO->GetIOComponentType();
>>>   const std::type_info &typeinfo typeInfo =
>>>imageIO->GetComponentTypeInfo();
>>>   // NumberOfComponents is usually one, but for non-scalar pixel
>>>types, it can be anything
>>>  unsigned int NumberOfComponents = imageIO->GetNumberOfComponents();
>>>  return 0;
>>>}
>>>
>>>Gaetan Lehmann wrote:
>>>      
>>>
>>>>Hi,
>>>>
>>>>Is it possible to write code able to read an image with an unknown
>>>>number of channels at build time, and to extract each channel in an
>>>>individual image ?
>>>>If yes, with which classes can it be done ?
>>>>
>>>>Thanks,
>>>>
>>>>Gaetan
>>>>
>>>>
>>>>------------------------------------------------------------------------
>>>>
>>>>_______________________________________________
>>>>Insight-users mailing list
>>>>Insight-users at itk.org
>>>>http://www.itk.org/mailman/listinfo/insight-users
>>>>        
>>>>
>>>_______________________________________________
>>>Insight-users mailing list
>>>Insight-users at itk.org
>>>http://www.itk.org/mailman/listinfo/insight-users
>>>      
>>>
>>>------------------------------------------------------------------------
>>>
>>>_______________________________________________
>>>Insight-users mailing list
>>>Insight-users at itk.org
>>>http://www.itk.org/mailman/listinfo/insight-users
>>>      
>>>


More information about the Insight-users mailing list