[Insight-users] separate interlaced frames

Luis Ibanez luis.ibanez at kitware.com
Mon Sep 3 18:29:18 EDT 2007


Hi Alex,

You shouldn't touch the spacing, origin, or orientation,
of images that are the output of a filter and are still
connected to the pipeline.

What you may want to use is the

        itkChangeImageInformationFilter
http://www.itk.org/Insight/Doxygen/html/classitk_1_1ChangeInformationImageFilter.html

By inserting this filter into your pipeline you will be
allowed you to change the spacing, origin, and orientation
of the image. Please make sure to provide call to both


    ImageType::SpacingType newSpacing = oldSpacing;
    newSpacing[1] *= 2;

    SetOutputSpacing( newSpacing )
    ChangeSpacingOn()



    Regards,


       Luis


---------------------------
Alex. GOUAILLARD wrote:
> hi all,
> 
> I am trying to separate two interlaced frames from a 1024*1024 image.
> 
> I modified the code from 
> Examples/Iterators/ImageLinearIteratorWithIndex.cxx
> to copy only one line out of two in the output image.
> It works.
> 
> But when I try to modify the spacing of the output image to be half the 
> one of the input, I have a memory access error thrown by the writer. I 
> must be missing something very simple, does anyone see it?
> 
> alex
> 
> //====================================================================
> #if defined(_MSC_VER)
> #pragma warning ( disable : 4786 )
> #endif
> 
> #include "itkImage.h"
> #include "itkRGBPixel.h"
> #include "itkImageLinearConstIteratorWithIndex.h"
> #include "itkImageLinearIteratorWithIndex.h"
> #include "itkImageFileReader.h"
> #include "itkImageFileWriter.h"
> 
> int main( int argc, char *argv[] )
> {
>   // Verify the number of parameters on the command line.
>   if ( argc < 3 )
>     {
>     std::cerr << "Missing parameters. " << std::endl;
>     std::cerr << "Usage: " << std::endl;
>     std::cerr << argv[0]
>               << " inputImageFile outputImageFile"
>               << std::endl;
>     return -1;
>     }
> 
>   // 2 to get one slice
>   // 3 to get on time point
>   // 4 if I get crazy enough one day to make one big file
>   const unsigned int Dimension = 2;
> 
>   typedef itk::RGBPixel< unsigned char > RGBPixelType;
>   typedef itk::Image< RGBPixelType, Dimension >  ImageType;
> 
>   typedef itk::ImageLinearIteratorWithIndex< ImageType > IteratorType;
>   typedef itk::ImageLinearConstIteratorWithIndex< ImageType > 
> ConstIteratorType;
> 
>   typedef itk::ImageFileReader< ImageType > ReaderType;
>   typedef itk::ImageFileWriter< ImageType > WriterType;
> 
>   ImageType::ConstPointer inputImage;
>   ReaderType::Pointer reader = ReaderType::New();
>   reader->SetFileName( argv[1] );
>   try
>     {
>     reader->Update();
>     inputImage = reader->GetOutput();
>     }
>   catch ( itk::ExceptionObject &err)
>     {
>     std::cout << "ExceptionObject caught a !" << std::endl;
>     std::cout << err << std::endl;
>     return -1;
>     }
> 
>   ImageType::Pointer outputImage1 = ImageType::New();
>   ImageType::RegionType Region =  inputImage->GetRequestedRegion();
>   ImageType::SizeType size = Region.GetSize( );
>   // size[1] = size[1] / 2; // only y: one line out of two
>   Region.SetSize( size );
> 
>   outputImage1->SetRegions( Region );
>   outputImage1->CopyInformation( inputImage );
>   outputImage1->Allocate();
>   ImageType::SpacingType spacing = inputImage->GetSpacing( );
>   //spacing[1] *= 2;
>   outputImage1->SetSpacing( spacing );
> 
>   ConstIteratorType inputIt( inputImage, 
> inputImage->GetRequestedRegion() );
>   IteratorType outputIt( outputImage1, inputImage->GetRequestedRegion() );
> 
>   inputIt.SetDirection(0);
>   outputIt.SetDirection(0);
> 
>   for ( inputIt.GoToBegin(),  outputIt.GoToBegin(); ! inputIt.IsAtEnd();
>         outputIt.NextLine(),  inputIt.NextLine())
>     {
>     inputIt.GoToBeginOfLine();
>     outputIt.GoToEndOfLine();
>     --outputIt;
>     while ( ! inputIt.IsAtEndOfLine() )
>       {
>       outputIt.Set( inputIt.Get() );
>       ++inputIt;
>       --outputIt;
>       }
>     // skip one line
>     inputIt.NextLine();
>     }
> 
>   WriterType::Pointer writer = WriterType::New();
>   writer->SetFileName( argv[2] );
>   writer->SetInput(outputImage1);
>   try
>     {
>     writer->Update();
>     }
>   catch ( itk::ExceptionObject &err)
>     {
>     std::cout << "ExceptionObject caught !" << std::endl;
>     std::cout << err << std::endl;
>     return -1;
>     }
> 
>   return 0;
> }
> 
> 
> 
> 
> _______________________________________________
> 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