[Insight-users] Dicom Series Read Series Write: change in intensity values

"Eliana M. Vásquez Osorio" e.vasquezosorio at erasmusmc.nl
Mon Nov 26 05:48:29 EST 2007


Hello,
I upgraded ITK and the problem is now solved.  The problem was related 
to the rescale slope and intersect tags in the DICOM files.

Thank you Mathieu for your help! ,
Eliana


Mathieu Malaterre wrote:
> That's correct the patch was added ~August 2007, so only ITK 3.4 and
> ITK 3.6 will have it.
>
> Sorry for the troubles,
> -Mathieu
>
> On Nov 22, 2007 3:42 PM, "Eliana M. Vásquez Osorio"
> <e.vasquezosorio at erasmusmc.nl> wrote:
>   
>> Hello Mathieu,
>> no.  I don't have it. -- I checked the whole file.  :\
>> Just checking the version of ITK I got installed... it's ITK 3.2...
>> :]  ups.
>> I think I'll upgrade the ITK version before I continue working...
>>
>>
>> Thank you,
>> Eliana
>>
>>
>> Mathieu Malaterre wrote:
>>     
>>> Hi Eliana,
>>>
>>>   Could you please check that you have something like this in your
>>> itkGDCMImageIO.cxx: (~line 835)
>>>
>>>           && !(dictEntry->GetGroup() == 0x0028 &&
>>> dictEntry->GetElement() == 0x1052)
>>>           && !(dictEntry->GetGroup() == 0x0028 &&
>>> dictEntry->GetElement() == 0x1053)
>>>
>>>   This means that itk-gdcm will *never* write a rescale slope/rescale
>>> intercept because of the issue you just mentionned.
>>>
>>> Thanks again for your feedback,
>>> -Mathieu
>>>
>>>
>>> On Nov 22, 2007 3:19 PM, "Eliana M. Vásquez Osorio"
>>> <e.vasquezosorio at erasmusmc.nl> wrote:
>>>
>>>       
>>>> [off list]
>>>>
>>>> Hello Mathieu,
>>>>
>>>> I made the change in the pixel type with no difference.  :(
>>>> But -- here comes the nice part --, while reading the last FAQ entry you
>>>> send, I became aware of the rescale slope and intercept.  And I know
>>>> what is happening now.
>>>>
>>>> I'm using the metadata of the read files into the written files (line
>>>> 92).  The rescale interpect is -1000 in both, the read and the written
>>>> image.  And that -1000 is exactly the difference I found between the
>>>> intensity values of the images (write - read).  So what happens is this:
>>>> - when reading the pixel value is calculated with [RWV = slope * SP +
>>>> intercept]  from FAQ
>>>> - when writing, the pixel value is not changed, and the rescale
>>>> intercept is set to -1000.
>>>> * when the written image is read, the pixel values are again calculated
>>>> and you get the intercept (-1000) added again.
>>>>
>>>> I read and write the same image 3 times, and the difference was accumulating
>>>> [original image] -> readNwrite -> [written1] -> readNwrite -> [written2]
>>>> -> readNwrite -> [written3]
>>>> and
>>>> written1 was -1000 offset from the original (pixel intensity values);
>>>> written2 was -2000 offset from the original (pixel intensity values);
>>>> written3 was -3000 offset from the original (pixel intensity values);
>>>>
>>>> I just include this in the code, and the problems seems to be fixed:
>>>>     WriterType::DictionaryArrayRawPointer dicts =
>>>> dcmreader->GetMetaDataDictionaryArray();
>>>>     WriterType::DictionaryArrayType::const_iterator itdict;
>>>>     for( itdict = dicts->begin(); itdict != dicts->end(); itdict++ )
>>>>     {
>>>>         WriterType::DictionaryRawPointer dicti = *itdict;
>>>>         //change 0028,1052,Rescale Intercept from -1000 to 0
>>>>         itk::EncapsulateMetaData<std::string>( *dicti, "0028|1052", "0" );
>>>>     }
>>>>
>>>> But the question that remains is, should this rescale intercept be
>>>> updated automatically when writing?  Are they other values that should
>>>> also be "corrected" when writing? Is this a bug? or is just a "special
>>>> case"?
>>>>
>>>> Thank you very much.  I'll write something to the list, just to round
>>>> off the 'topic'.
>>>> (sorry for so much ping pong)
>>>>
>>>> Eliana
>>>>
>>>>
>>>> Mathieu Malaterre wrote:
>>>>
>>>>         
>>>>> Hi Eliana,
>>>>>
>>>>>   This is common mistake, I even do it myself. ITK toolkit requires
>>>>> you to instanciate using the proper pixel type. In this case you
>>>>> should have used signed short and not unsigned short.
>>>>>   I have added a FAQ entry:
>>>>>
>>>>> http://www.itk.org/Wiki/ITK_FAQ#DICOM:_Rescale_Slope.2FIntercept
>>>>>
>>>>>   I have no other better solution for now than explanation. At some
>>>>> point I might trigger an exception in the reader to indicate that the
>>>>> operation might be invalid. But I have to see if ITK-dev are ok with
>>>>> that (this will break ITK design I think).
>>>>>
>>>>> -Mathieu
>>>>>
>>>>> On Nov 22, 2007 11:33 AM, "Eliana M. Vásquez Osorio"
>>>>> <e.vasquezosorio at erasmusmc.nl> wrote:
>>>>>
>>>>>
>>>>>           
>>>>>> [off list]
>>>>>>
>>>>>> Hello Mathieu!
>>>>>> Thank you very much! :) :) :)
>>>>>> I attach the code (.cxx), two "to-be-read" images, and two written
>>>>>> images (in zip files).
>>>>>> The corresponding files are:
>>>>>> ct344.dcm (read.zip) == ct000.dcm (written.zip)
>>>>>> ct345.dcm (read.zip) == ct001.dcm (written.zip)
>>>>>>
>>>>>> Have a nice day,
>>>>>> Eliana
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>> Mathieu Malaterre wrote:
>>>>>>
>>>>>>
>>>>>>             
>>>>>>> [off list]
>>>>>>>
>>>>>>> Hi Eliana,
>>>>>>>
>>>>>>>   Yes could you please send me your c++ code (also I trust you when
>>>>>>> you say you are not applying any filter). And at the same time, send
>>>>>>> me at least two files, the one read and the equivalent when written.
>>>>>>> It sill save me time to look into your issue.
>>>>>>>
>>>>>>> Thanks
>>>>>>> -Mathieu
>>>>>>>
>>>>>>>
>>>>>>>
>>>>>>> On Nov 21, 2007 3:15 PM, "Eliana M. Vásquez Osorio"
>>>>>>> <e.vasquezosorio at erasmusmc.nl> wrote:
>>>>>>>
>>>>>>>
>>>>>>>
>>>>>>>               
>>>>>>>> Hello Mathieu,
>>>>>>>> I have read it and it makes sense  :), but I think it does not help me yet.
>>>>>>>>  From what I read, the written image should keep the intensity values
>>>>>>>> equal. (since the XXX from bit 15 to bit 12 are set to 0).  I do not
>>>>>>>> have to apply any filter... right?
>>>>>>>> But, still the intensities in the images are different.  :(
>>>>>>>> Should I include the code or upload images to make things clearer?
>>>>>>>>
>>>>>>>> Thank you,
>>>>>>>> Eliana
>>>>>>>>
>>>>>>>>
>>>>>>>> Mathieu Malaterre wrote:
>>>>>>>>
>>>>>>>>
>>>>>>>>
>>>>>>>>                 
>>>>>>>>> On Nov 20, 2007 6:13 PM, "Eliana M. Vásquez Osorio"
>>>>>>>>> <e.vasquezosorio at erasmusmc.nl> wrote:
>>>>>>>>>
>>>>>>>>>
>>>>>>>>>
>>>>>>>>>
>>>>>>>>>                   
>>>>>>>>>> Hello itk-users,
>>>>>>>>>> I have a problem when writing (?) images series.  The intensity values
>>>>>>>>>> change.
>>>>>>>>>> This is what I did:
>>>>>>>>>> - Read a series of dicom image into a volume
>>>>>>>>>> - Write the volume out as another series of dicom files in a different
>>>>>>>>>> folder.
>>>>>>>>>> Since I do nothing in between (no filter is applied), I expected to get
>>>>>>>>>> the same images after the writing process, but I found differences in
>>>>>>>>>> the intensity values between the images.
>>>>>>>>>>
>>>>>>>>>> My input images use 12 bits for storing image info. (Bits Stored=12,
>>>>>>>>>> 0028|0101 dicom tag).  The written image uses 16 bits.  I imagine that
>>>>>>>>>> the intensity changes are related to this... but I don't get a clear
>>>>>>>>>> idea how/why this changes the image that is produced.  Looking closely
>>>>>>>>>> to the data, I found an offset of 1000 in the written images (compared
>>>>>>>>>> to the original images).  Should I apply a filter to "correct" the 12
>>>>>>>>>> bit to a 16 bit?  if so, which filter?  a casting filter?
>>>>>>>>>>
>>>>>>>>>> Any idea or suggestion is welcome! :)
>>>>>>>>>>
>>>>>>>>>>
>>>>>>>>>>
>>>>>>>>>>
>>>>>>>>>>                     
>>>>>>>>> This is in the  FAQ:
>>>>>>>>>
>>>>>>>>> http://www.itk.org/Wiki/ITK_FAQ#DICOM:_Bits_Allocated.2C_Bits_Stored_and_High_Bit
>>>>>>>>>
>>>>>>>>> Let me know if this makes any sense :)
>>>>>>>>>
>>>>>>>>>
>>>>>>>>>
>>>>>>>>>
>>>>>>>>>
>>>>>>>>>                   
>>>>>>>               
>>>>>> --
>>>>>> Eliana M. Vásquez Osorio
>>>>>> PhD Student (OIO)
>>>>>> Eramus MC - Daniel den Hoed
>>>>>> e.vasquezosorio at erasmusmc.nl
>>>>>>
>>>>>> Radiotherapie, G2-75
>>>>>> Erasmus MC - Daniel den Hoed
>>>>>> Groene Hilledijk 301,
>>>>>> 3075 EA Rotterdam
>>>>>> T: +31 10 – 43 91 491
>>>>>> F: +31 10 – 43 91 012
>>>>>>
>>>>>>
>>>>>> #if defined(_MSC_VER)
>>>>>> #pragma warning ( disable : 4786 )
>>>>>> #endif
>>>>>>
>>>>>> #ifdef __BORLANDC__
>>>>>> #define ITK_LEAN_AND_MEAN
>>>>>> #endif
>>>>>>
>>>>>>
>>>>>> #include "itkGDCMImageIO.h"
>>>>>> #include "itkGDCMSeriesFileNames.h"
>>>>>> #include "itkImageSeriesReader.h"
>>>>>> #include "itkWarpImageFilter.h"
>>>>>> #include "itkLinearInterpolateImageFunction.h"
>>>>>> #include "itkImageSeriesWriter.h"
>>>>>> #include "itkNumericSeriesFileNames.h"
>>>>>> #include "itkShiftScaleImageFilter.h"
>>>>>>
>>>>>> int main( int argc, char* argv[] )
>>>>>> {
>>>>>>         if( argc < 4 )
>>>>>>         {
>>>>>>                 std::cerr << "Usage: " << std::endl;
>>>>>>                 std::cerr << argv[0] << " DicomDirectory  vectorFieldAsImageFile FolderForNewImages"
>>>>>>                         << std::endl;
>>>>>>                 return EXIT_FAILURE;
>>>>>>         }
>>>>>>
>>>>>>         // Read the dicom series
>>>>>>         typedef unsigned short    PixelType;
>>>>>>         const unsigned int      Dimension = 3;
>>>>>>         typedef itk::Image<PixelType, Dimension>  DCMImageType;
>>>>>>         typedef itk::ImageSeriesReader<DCMImageType> DCMReaderType;
>>>>>>         typedef itk::GDCMImageIO ImageIOType;
>>>>>>
>>>>>>         DCMReaderType::Pointer dcmreader = DCMReaderType::New();
>>>>>>         ImageIOType::Pointer dicomIO = ImageIOType::New();
>>>>>>         dcmreader->SetImageIO( dicomIO );
>>>>>>
>>>>>>         typedef itk::GDCMSeriesFileNames DCMNamesGeneratorType;
>>>>>>         DCMNamesGeneratorType::Pointer nameGenerator = DCMNamesGeneratorType::New();
>>>>>>         nameGenerator->SetUseSeriesDetails( true );
>>>>>>         nameGenerator->SetDirectory( argv[1] );
>>>>>>
>>>>>>         try
>>>>>>         {
>>>>>>                 typedef std::vector< std::string >    SeriesIdContainer;
>>>>>>                 const SeriesIdContainer & seriesUID = nameGenerator->GetSeriesUIDs();
>>>>>>                 if( seriesUID.size() > 1 )
>>>>>>                         std::cout << std::endl << "The directory: " << std::endl << std::endl << argv[1] << std::endl << std::endl
>>>>>>                                           << "Contains " << seriesUID.size() <<  " DICOM Series. Using 1st. serie" << std::endl << std::endl;
>>>>>>
>>>>>>                 std::string seriesIdentifier = seriesUID.begin()->c_str();
>>>>>>                 std::cout << std::endl << std::endl << "Now reading series: " << std::endl << std::endl
>>>>>>                               << seriesIdentifier << std::endl << std::endl << std::endl;
>>>>>>
>>>>>>                 typedef std::vector< std::string >   FileNamesContainer;
>>>>>>                 FileNamesContainer fileNames;
>>>>>>                 fileNames = nameGenerator->GetFileNames( seriesIdentifier );
>>>>>>
>>>>>>                 dcmreader->SetFileNames( fileNames );
>>>>>>                 try{ dcmreader->Update(); }
>>>>>>                 catch (itk::ExceptionObject &ex) {      std::cout << ex << std::endl;
>>>>>> return EXIT_FAILURE; }
>>>>>>         }
>>>>>>         catch (itk::ExceptionObject &ex) { std::cout << ex << std::endl;
>>>>>> return EXIT_FAILURE; }
>>>>>>
>>>>>>         typedef itk::Image<PixelType, 2> Image2DType;
>>>>>>         typedef itk::ImageSeriesWriter< DCMImageType, Image2DType > WriterType;
>>>>>>         WriterType::Pointer writer = WriterType::New();
>>>>>>         writer->SetInput( dcmreader->GetOutput() );
>>>>>>
>>>>>>         typedef itk::NumericSeriesFileNames    NumNameGeneratorType;
>>>>>>         NumNameGeneratorType::Pointer numNameGenerator = NumNameGeneratorType::New();
>>>>>>         std::string format = std::string(argv[3]) + "/ct%03d.dcm";   // filename extension
>>>>>>         numNameGenerator->SetSeriesFormat( format.c_str() );
>>>>>>         DCMImageType::ConstPointer inputImage = dcmreader->GetOutput();
>>>>>>         DCMImageType::RegionType   region     = inputImage->GetLargestPossibleRegion();
>>>>>>         DCMImageType::IndexType    start      = region.GetIndex();
>>>>>>         DCMImageType::SizeType     size       = region.GetSize();
>>>>>>
>>>>>>         const unsigned int firstSlice = start[2];
>>>>>>         const unsigned int lastSlice  = start[2] + size[2] - 1;
>>>>>>
>>>>>>         numNameGenerator->SetStartIndex( firstSlice );
>>>>>>         numNameGenerator->SetEndIndex( lastSlice );
>>>>>>         numNameGenerator->SetIncrementIndex( 1 );
>>>>>>         writer->SetFileNames( numNameGenerator->GetFileNames() );
>>>>>>
>>>>>>         writer->SetImageIO( dicomIO );
>>>>>>         writer->SetMetaDataDictionaryArray( dcmreader->GetMetaDataDictionaryArray() );
>>>>>>
>>>>>>         try { writer->Update(); }
>>>>>>         catch (itk::ExceptionObject &ex) { std::cout << ex << std::endl;
>>>>>> return EXIT_FAILURE; }
>>>>>>
>>>>>>         return EXIT_SUCCESS;
>>>>>> }
>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>>             
>>>>>
>>>>>           
>>>> --
>>>> Eliana M. Vásquez Osorio
>>>> PhD Student (OIO)
>>>> Eramus MC - Daniel den Hoed
>>>> e.vasquezosorio at erasmusmc.nl
>>>>
>>>> Radiotherapie, G2-75
>>>> Erasmus MC - Daniel den Hoed
>>>> Groene Hilledijk 301,
>>>> 3075 EA Rotterdam
>>>> T: +31 10 – 43 91 491
>>>> F: +31 10 – 43 91 012
>>>>
>>>>
>>>>
>>>>         
>>>
>>>
>>>       
>> --
>> Eliana M. Vásquez Osorio
>> PhD Student (OIO)
>> Eramus MC - Daniel den Hoed
>> e.vasquezosorio at erasmusmc.nl
>>
>> Radiotherapie, G2-75
>> Erasmus MC - Daniel den Hoed
>> Groene Hilledijk 301,
>> 3075 EA Rotterdam
>> T: +31 10 – 43 91 491
>> F: +31 10 – 43 91 012
>>
>>
>>     
>
>
>
>   


-- 
Eliana M. Vásquez Osorio
PhD Student (OIO)
Eramus MC - Daniel den Hoed
e.vasquezosorio at erasmusmc.nl

Radiotherapie, G2-75
Erasmus MC - Daniel den Hoed
Groene Hilledijk 301, 
3075 EA Rotterdam
T: +31 10 – 43 91 491 
F: +31 10 – 43 91 012



More information about the Insight-users mailing list