[Insight-users] Simple resampling problem

Richard Beare richard.beare at gmail.com
Thu Feb 18 18:33:54 EST 2010


OK,
Here's a series of resampling tests:
Test 1 produces the expected result, but tests 2 and 3 produce blank images.

Test 1) Upsampling the original image:


Input image:

Image (0x2f72570)
  RTTI typeinfo:   itk::Image<short, 3u>
  Reference Count: 3
  Modified Time: 189
  Debug: Off
  Observers:
    none
  Source: (none)
  Source output index:  0
  Release Data: Off
  Data Released: False
  Global Release Data: Off
  PipelineMTime: 0
  UpdateMTime: 182
  LargestPossibleRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [256, 176, 20]
  BufferedRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [256, 176, 20]
  RequestedRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [256, 176, 20]
  Spacing: [1.25, 1.25, 6]
  Origin: [-165.833, 127.35, 80.3368]
  Direction:
0.989374 0 0.145394
0.026043 -0.983827 -0.177216
-0.143043 -0.17912 0.973373

  IndexToPointMatrix:
  1.23672 0 0.872366
0.0325538 -1.22978 -1.0633
-0.178804 -0.2239 5.84024

  PointToIndexMatrix:
  0.791499 0.0208344 -0.114434
1.17893e-08 -0.787062 -0.143296
0.0242324 -0.0295361 0.162229

  PixelContainer:
    ImportImageContainer (0x2f731b0)
      RTTI typeinfo:   itk::ImportImageContainer<unsigned long, short>
      Reference Count: 1
      Modified Time: 179
      Debug: Off
      Observers:
        none
      Pointer: 0x3800000
      Container manages memory: true
      Size: 901120
      Capacity: 901120

Output Image:

Image (0x2f74750)
  RTTI typeinfo:   itk::Image<short, 3u>
  Reference Count: 2
  Modified Time: 245
  Debug: Off
  Observers:
    none
  Source: (none)
  Source output index:  0
  Release Data: Off
  Data Released: False
  Global Release Data: Off
  PipelineMTime: 0
  UpdateMTime: 238
  LargestPossibleRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [256, 176, 40]
  BufferedRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [256, 176, 40]
  RequestedRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [256, 176, 40]
  Spacing: [1.25, 1.25, 3]
  Origin: [-165.833, 127.35, 80.3368]
  Direction:
0.989374 0 0.145394
0.026043 -0.983827 -0.177216
-0.143043 -0.17912 0.973373

  IndexToPointMatrix:
  1.23672 0 0.436183
0.0325538 -1.22978 -0.531649
-0.178804 -0.2239 2.92012

  PointToIndexMatrix:
  0.791499 0.0208344 -0.114434
1.17893e-08 -0.787062 -0.143296
0.0484648 -0.0590721 0.324458

  PixelContainer:
    ImportImageContainer (0x2f74290)
      RTTI typeinfo:   itk::ImportImageContainer<unsigned long, short>
      Reference Count: 1
      Modified Time: 236
      Debug: Off
      Observers:
        none
      Pointer: 0x39b8000
      Container manages memory: true
      Size: 1802240
      Capacity: 1802240


Test 2) Upsampling of the extracted image

Input Image:Image (0x2f75540)
  RTTI typeinfo:   itk::Image<short, 3u>
  Reference Count: 2
  Modified Time: 220
  Debug: Off
  Observers:
    none
  Source: (none)
  Source output index:  0
  Release Data: Off
  Data Released: False
  Global Release Data: Off
  PipelineMTime: 0
  UpdateMTime: 213
  LargestPossibleRegion:
    Dimension: 3
    Index: [112, 93, 0]
    Size: [40, 40, 20]
  BufferedRegion:
    Dimension: 3
    Index: [112, 93, 0]
    Size: [40, 40, 20]
  RequestedRegion:
    Dimension: 3
    Index: [112, 93, 0]
    Size: [40, 40, 20]
  Spacing: [1.25, 1.25, 6]
  Origin: [-165.833, 127.35, 80.3368]
  Direction:
0.989374 0 0.145394
0.026043 -0.983827 -0.177216
-0.143043 -0.17912 0.973373

  IndexToPointMatrix:
  1.23672 0 0.872366
0.0325538 -1.22978 -1.0633
-0.178804 -0.2239 5.84024

  PointToIndexMatrix:
  0.791499 0.0208344 -0.114434
1.17893e-08 -0.787062 -0.143296
0.0242324 -0.0295361 0.162229

  PixelContainer:
    ImportImageContainer (0x2f74330)
      RTTI typeinfo:   itk::ImportImageContainer<unsigned long, short>
      Reference Count: 1
      Modified Time: 211
      Debug: Off
      Observers:
        none
      Pointer: 0x2eaf000
      Container manages memory: true
      Size: 32000
      Capacity: 32000

Output Image:

Image (0x2f74890)
  RTTI typeinfo:   itk::Image<short, 3u>
  Reference Count: 2
  Modified Time: 264
  Debug: Off
  Observers:
    none
  Source: (none)
  Source output index:  0
  Release Data: Off
  Data Released: False
  Global Release Data: Off
  PipelineMTime: 0
  UpdateMTime: 255
  LargestPossibleRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [40, 40, 40]
  BufferedRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [40, 40, 40]
  RequestedRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [40, 40, 40]
  Spacing: [1.25, 1.25, 3]
  Origin: [-165.833, 127.35, 80.3368]
  Direction:
0.989374 0 0.145394
0.026043 -0.983827 -0.177216
-0.143043 -0.17912 0.973373

  IndexToPointMatrix:
  1.23672 0 0.436183
0.0325538 -1.22978 -0.531649
-0.178804 -0.2239 2.92012

  PointToIndexMatrix:
  0.791499 0.0208344 -0.114434
1.17893e-08 -0.787062 -0.143296
0.0484648 -0.0590721 0.324458

  PixelContainer:
    ImportImageContainer (0x2f74a60)
      RTTI typeinfo:   itk::ImportImageContainer<unsigned long, short>
      Reference Count: 1
      Modified Time: 253
      Debug: Off
      Observers:
        none
      Pointer: 0x2ebf000
      Container manages memory: true
      Size: 64000
      Capacity: 64000

Test 3) Feed extracted image through the same steps but leave spacing the same

Image (0x2f75540)
  RTTI typeinfo:   itk::Image<short, 3u>
  Reference Count: 2
  Modified Time: 220
  Debug: Off
  Observers:
    none
  Source: (none)
  Source output index:  0
  Release Data: Off
  Data Released: False
  Global Release Data: Off
  PipelineMTime: 0
  UpdateMTime: 213
  LargestPossibleRegion:
    Dimension: 3
    Index: [112, 93, 0]
    Size: [40, 40, 20]
  BufferedRegion:
    Dimension: 3
    Index: [112, 93, 0]
    Size: [40, 40, 20]
  RequestedRegion:
    Dimension: 3
    Index: [112, 93, 0]
    Size: [40, 40, 20]
  Spacing: [1.25, 1.25, 6]
  Origin: [-165.833, 127.35, 80.3368]
  Direction:
0.989374 0 0.145394
0.026043 -0.983827 -0.177216
-0.143043 -0.17912 0.973373

  IndexToPointMatrix:
  1.23672 0 0.872366
0.0325538 -1.22978 -1.0633
-0.178804 -0.2239 5.84024

  PointToIndexMatrix:
  0.791499 0.0208344 -0.114434
1.17893e-08 -0.787062 -0.143296
0.0242324 -0.0295361 0.162229

  PixelContainer:
    ImportImageContainer (0x2f74330)
      RTTI typeinfo:   itk::ImportImageContainer<unsigned long, short>
      Reference Count: 1
      Modified Time: 211
      Debug: Off
      Observers:
        none
      Pointer: 0x2eaf000
      Container manages memory: true
      Size: 32000
      Capacity: 32000

Image (0x2f74890)
  RTTI typeinfo:   itk::Image<short, 3u>
  Reference Count: 2
  Modified Time: 262
  Debug: Off
  Observers:
    none
  Source: (none)
  Source output index:  0
  Release Data: Off
  Data Released: False
  Global Release Data: Off
  PipelineMTime: 0
  UpdateMTime: 253
  LargestPossibleRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [40, 40, 20]
  BufferedRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [40, 40, 20]
  RequestedRegion:
    Dimension: 3
    Index: [0, 0, 0]
    Size: [40, 40, 20]
  Spacing: [1.25, 1.25, 6]
  Origin: [-165.833, 127.35, 80.3368]
  Direction:
0.989374 0 0.145394
0.026043 -0.983827 -0.177216
-0.143043 -0.17912 0.973373

  IndexToPointMatrix:
  1.23672 0 0.872366
0.0325538 -1.22978 -1.0633
-0.178804 -0.2239 5.84024

  PointToIndexMatrix:
  0.791499 0.0208344 -0.114434
1.17893e-08 -0.787062 -0.143296
0.0242324 -0.0295361 0.162229

  PixelContainer:
    ImportImageContainer (0x2f74a60)
      RTTI typeinfo:   itk::ImportImageContainer<unsigned long, short>
      Reference Count: 1
      Modified Time: 251
      Debug: Off
      Observers:
        none
      Pointer: 0x2ebf000
      Container manages memory: true
      Size: 32000
      Capacity: 32000


On Wed, Feb 17, 2010 at 2:04 AM, Luis Ibanez <luis.ibanez at kitware.com> wrote:
> Hi Richard,
>
> If you set
>
>          NewSpacing = inputSpacing,
>
> Do you get as output an image
> that looks like the input image ?
>
>
> What are the values of
>
> * inputSpacing ?
> * inputSize ?
> * inputOrigin ?
> * inputDirection ?
>
>
> Please let us know,
>
>
>    Thanks
>
>
>         Luis
>
>
> ------------------------------------------
> On Tue, Feb 16, 2010 at 1:21 AM, Richard Beare <richard.beare at gmail.com> wrote:
>> Hi,
>> I'm using an obviously buggy resampling procedure to upsample a small
>> image that has been extracted from a larger one, and hence has non
>> zero origin information. My result is blank, so I'm missing a setting
>> somewhere, but don't have a clue as to what. Can anyone spot the
>> problem?
>>
>> template  <class RawIm>
>> typename RawIm::Pointer upsampleIm(typename RawIm::Pointer input,
>> typename RawIm::SpacingType NewSpacing, int interp=1)
>> {
>>  const int dim = RawIm::ImageDimension;
>>  typedef typename RawIm::PixelType PixelType;
>>
>>  typedef typename itk::ResampleImageFilter<RawIm, RawIm >  ResampleFilterType;
>>  typedef typename itk::IdentityTransform< double, dim >  TransformType;
>>  typename ResampleFilterType::Pointer resampler = ResampleFilterType::New();
>>
>>  input->Update();
>>
>>  typename TransformType::Pointer transform = TransformType::New();
>>  transform->SetIdentity();
>>  resampler->SetTransform( transform );
>>  typedef typename itk::LinearInterpolateImageFunction<RawIm, double >
>>  LInterpolatorType;
>>  typedef typename itk::NearestNeighborInterpolateImageFunction<RawIm,
>> double >  NNInterpolatorType;
>>
>>  typename ResampleFilterType::InterpolatorPointerType interpolator;
>>  switch (interp)
>>    {
>>    case 0:
>>      interpolator = NNInterpolatorType::New();
>>      break;
>>    case 1:
>>      interpolator = LInterpolatorType::New();
>>      break;
>>    default:
>>      std::cout << "Unsupported interpolator" << std::endl;
>>    }
>>
>>  resampler->SetInterpolator( interpolator );
>>  resampler->SetDefaultPixelValue( 0 );
>>
>>  const typename RawIm::SpacingType& inputSpacing = input->GetSpacing();
>>  typename RawIm::SpacingType spacing;
>>  typename RawIm::SizeType   inputSize =
>> input->GetLargestPossibleRegion().GetSize();
>>  typename RawIm::SizeType   size;
>>  typedef typename RawIm::SizeType::SizeValueType SizeValueType;
>>
>>
>>  for (int i = 0; i < dim; i++)
>>    {
>>    float factor = inputSpacing[i]/NewSpacing[i];
>>    size[i] = static_cast< SizeValueType >( inputSize[i] * factor );
>>    }
>>   std::cout << inputSpacing << NewSpacing << std::endl;
>>   std::cout << inputSize << size << input->GetOrigin() << std::endl;
>>
>>  resampler->SetSize( size );
>>  resampler->SetOutputSpacing( NewSpacing );
>>  resampler->SetOutputOrigin( input->GetOrigin() );
>>  resampler->SetOutputDirection(input->GetDirection());
>>  resampler->SetInput(input);
>>  typename RawIm::Pointer result = resampler->GetOutput();
>>  result->Update();
>>  result->DisconnectPipeline();
>>  return(result);
>> }
>>
>> //////////////////////////////////////////////////
>> _____________________________________
>> Powered by www.kitware.com
>>
>> Visit other Kitware open-source projects at
>> http://www.kitware.com/opensource/opensource.html
>>
>> Kitware offers ITK Training Courses, for more information visit:
>> http://www.kitware.com/products/protraining.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
>>
>


More information about the Insight-users mailing list