[Insight-users] Nifti vs Dicom orientation

Luke Bloy luke.bloy at gmail.com
Mon Apr 6 12:13:20 EDT 2009


Hi all,

So this is my understanding of how orientations and ITK work.

the GetDirection and GetOrigin methods in itk::Image return the 
directional cosines in an LPS world frame. Thus an identity matrix 
returned by GetDirection() would mean the image is stored in LPS.

The QForm and SForm matrices in the nifti file format  are in the RAS 
world coordinate frame. See 
http://nifti.nimh.nih.gov/nifti-1/documentation/nifti1fields/nifti1fields_pages/qsform.html

So Mathieu, what you are reporting is the expected behavior.
qform (shown by fslhd etc) of
[1, 0, 0, -90]
[0, 1, 0, -126]
[0, 0, 1, -72]
[0, 0, 0, 1]

should yield a GetDirection() and GetOrigin of
-1  0  0   90
0  -1  0   126
0   0  1  -72
0   0  0   1


Antoine,
Your situation is more confusing, I'm not familiar with dcm2nii, but it 
isn't unheard of to have dicom converters change the orientation of the 
output image, to bring it into some consistent frame (LAS is a favorite 
of labs that have lots of old code expecting analyze files). From what 
you posted it seems that dcm2nii converted your image from an LPS 
arranged dicom to an LAS arranged Nifti file. I would investigate the 
nifti file's qform and sform matrics using either nifti_tool, fslhd or 
pynifti to confirm this.

Hope this helps,
Luke

Mathieu Coursolle wrote:
> Hi,
>
> I am having a similar issue, which may be related.
>
> I have a NIfTI file, which I load using itkNIfTIImageIO.
>
> By stepping into ITK, I can see that the s-form (sto_xyz) is:
>
> [1, 0, 0, -90]
> [0, 1, 0, -126]
> [0, 0, 1, -72]
> [0, 0, 0, 1]
>
> However, when it read the origin and orientation of the resulting 
> itkImage,
> I get:
>
> [-1, 0, 0, 90]
> [0, -1, 0, 126]
> [0, 0, 1, -72]
> [0, 0, 0, 1]
>
> This data is in the MNI-152 space. I know that the index (0, 0, 0) 
> should be
> (-90, -126, -72) mm in LPI  coordinate system (NIfTI's default I believe).
>
> Since I want my data in the RAI coordinate system, I then use the 
> itkOrientImageFilter.
> Assuming it converts the image from LPI to RAI, I would expect the 
> ouput to be:
> [-1, 0, 0, 90]
> [0, -1, 0, 90]
> [0, 0, 1, -72]
> [0, 0, 0, 1]
>
> But I get:
> [1, 0, 0, -90]
> [0, 1, 0, -90]
> [0, 0, 1, -72]
> [0, 0, 0, 1]
>
> By the way, my image dimensions are 181 x 217 x 181.
>
> I know that index (0, 0, 0) is (90,  90, -72) mm in RAI. 
>
> I am right that the matrix I get from ITK is not what is expected ?
>
> It looks to me that the itkNIfTIImageIO is doing some transform to
> the orientation matrix.
>
> What would be that transformation's purpose ?
>
> Thanks a lot!
>
> Mathieu
>
> ____________________ 
> Mathieu Coursolle, M.Ing.
> Rogue Research Inc.        
> www.rogue-research.com <http://www.rogue-research.com/>
>
>
> On 6-Apr-09, at 10:46 AM, Antoine DUCHAMPS wrote:
>
>> Hi all,
>>
>> I have a DWI sequence (brain) acquired with a Siemens scanner in Dicom
>> format. After converting it to NIfTI with dcm2nii
>> (http://www.sph.sc.edu/comd/rorden/mricron/dcm2nii.html) I tried to read
>> it an recover the orientation by using ITK (I've copied the code below).
>> The matrix I obtain is the following:
>>
>> 0.9997  0.0000  0.0252
>> 0.0002 -1.0000 -0.0072
>> -0.0252 -0.0072  0.9997
>>
>> However, the transforation matrix in the Dicom header is
>>
>> 0.9997  0.0000  0.0252
>> 0.0002  1.0000 -0.0072
>> -0.0252  0.0072  0.9997
>>
>> If the Dicom image orientation is LPS, this means that the image
>> orientation in NIfTI is LAS (And not RAS as I believed). Could anybody
>> clarify this please? Are there several possible orientations in NIfTI?
>> If so, how can I know the specific image orientation?
>>
>> Antoine.  
>>
>>
>> #include "itkImageFileReader.h"
>> #include "itkOrientedImage.h"
>>
>> const unsigned int Dimension = 4;
>> typedef short PixelType;
>> typedef itk::OrientedImage<PixelType, Dimension> ImageType;
>> typedef ImageType::Pointer ImagePointerType;
>> typedef itk::ImageFileReader< ImageType > ReaderType;
>>
>>
>> int main(int argc, char* argv[])
>> {
>>
>>  ImageType::Pointer image = ImageType::New();
>>
>>  ReaderType::Pointer reader = ReaderType::New();
>>  reader->SetFileName(argv[1]);
>>  reader->Update();
>>
>>  image = reader -> GetOutput();
>>  std::cout << image->GetDirection() << std::endl; std::cout.flush();
>>
>> }
>>
>> _____________________________________
>> Powered by www.kitware.com <http://www.kitware.com>
>>
>> Visit other Kitware open-source projects at
>> http://www.kitware.com/opensource/opensource.html
>>
>> Please keep messages on-topic and check the ITK FAQ at: 
>> http://www.itk.org/Wiki/ITK_FAQ
>>
>> Follow this link to subscribe/unsubscribe:
>> http://www.itk.org/mailman/listinfo/insight-users
>>
>
>
>
>
>
>
>
>
> ------------------------------------------------------------------------
>
> _____________________________________
> Powered by www.kitware.com
>
> Visit other Kitware open-source projects at
> http://www.kitware.com/opensource/opensource.html
>
> Please keep messages on-topic and check the ITK FAQ at: http://www.itk.org/Wiki/ITK_FAQ
>
> Follow this link to subscribe/unsubscribe:
> http://www.itk.org/mailman/listinfo/insight-users
>   
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://www.itk.org/pipermail/insight-users/attachments/20090406/36529814/attachment.htm>


More information about the Insight-users mailing list