View Issue Details Jump to Notes ] Print ]
IDProjectCategoryView StatusDate SubmittedLast Update
0012684VTK(No Category)public2011-10-28 12:022013-12-19 07:24
ReporterJohn Stark 
Assigned ToBill Lorensen 
PrioritynormalSeverityminorReproducibilityhave not tried
StatusclosedResolutionfixed 
PlatformOSOS Version
Product Version 
Target VersionFixed in Version6.1.0 
Summary0012684: vtkImageImport has problems with PipelineMtime and upstream changes.
DescriptionI have been using the ImageToVTKImageFilter, but have found a bug that the downstream pipeline does not correctly detect changes upstream of the filter.

In the sample below, an ITK image is constructed and filled with the value 127. This value is correctly extracted as float using a vtkImageCast filter. However, when I change the data buffer, filling it with zeros, the output initially reports the image value as unchanged, even after Update() is called. It appears to explicitly require the ITK converter to be updated, which should not be necessary.

I think the bug relates to vtkImageImport, which does not correctly report the PipelineMtime, but I am not sure.

I am using a git version of VTK, based on dbb88cda6b6d703b237b431629562ff275aae000 from July this year.

Code sample:

#include <itkImage.h>
#include <itkImageToVTKImageFilter.h>

#include <vtkImageCast.h>
#include <vtkImageData.h>
#include <vtkSmartPointer.h>

int main()
{
  const int dims[3] = { 64, 64, 64 }; // Image size in pixels
  
  typedef itk::Image<unsigned char, 3> ImageType;

  ImageType::RegionType region;
  for ( int i(0); i<3; ++i ) {
      region.SetSize(i, dims[i] );
      region.SetIndex(i,0);
  }

  ImageType::Pointer image = ImageType::New();
  image->SetLargestPossibleRegion(region);
  image->SetBufferedRegion(region);
  image->Allocate();
  image->FillBuffer( 127 );


  typedef itk::ImageToVTKImageFilter< ImageType > ConverterType;
  ConverterType::Pointer myConverter = ConverterType::New();
  myConverter->SetInput ( image );
  myConverter->UpdateOutputInformation();

  vtkSmartPointer<vtkImageCast> imCast = vtkSmartPointer<vtkImageCast>::New();
  imCast->SetOutputScalarTypeToFloat();
  imCast->SetInput(myConverter->GetOutput());

  vtkImageData * tstImg = imCast->GetOutput();
  tstImg->Update();
  float fval0 = ((const float*)tstImg->GetScalarPointer())[0];

  image->FillBuffer(0);
  image->Modified();
  image->GetPixelContainer()->Modified();

  tstImg->Update();
  float fval1 = ((const float*)tstImg->GetScalarPointer())[0];

  myConverter->GetOutput()->UpdateInformation();
  tstImg->Update();
  float fval2 = ((const float*)tstImg->GetScalarPointer())[0];

  if ( (fval0 == 127.)
      && (fval1 == 0.)
      && (fval2 == 0.) ) {
    return EXIT_SUCCESS;
  } else {
      return EXIT_FAILURE;
  }
}
TagsNo tags attached.
ProjectTBD
Typeincorrect functionality
Attached Filespatch file icon ComputePipelineMTime.patch [^] (1,746 bytes) 2012-11-06 08:23 [Show Content]
cxx file icon sampleCrash.cxx [^] (4,387 bytes) 2012-11-06 08:38
txt file icon CMakeLists.txt [^] (637 bytes) 2013-12-18 17:13 [Show Content]
cxx file icon main.cxx [^] (4,392 bytes) 2013-12-18 17:13

 Relationships

  Notes
(0027620)
John Stark (reporter)
2011-10-28 12:07

One possible solution may be to implement vtkAlgorithm::ComputePipelineMTime in vtkImageImport to return this->GetMTime(), but I have not tested this.
(0029655)
John Stark (reporter)
2012-11-06 08:31
edited on: 2012-12-12 07:52

I have attached a patch for vtk5.10.0 which indeed fixes this bug as suggested above.

This bug can lead to a crash if the input is changed to a new image, and the old image is deleted. This can leave the down stream pipeline pointing to stale memory, usually causing a seg fault.

A sample of how this can happen is also attached.

I pushed a patch to gerrit : http://review.source.kitware.com/#/c/8885/ [^]

(0029870)
John Stark (reporter)
2012-12-12 07:55
edited on: 2013-12-16 16:41

This also fixed problems for other users : http://www.vtk.org/pipermail/vtkusers/2012-November/126372.html [^]

Edit: Update link: http://public.kitware.com/pipermail/vtkusers/2012-November/077414.html [^]

(0032026)
John Stark (reporter)
2013-12-18 17:16
edited on: 2013-12-18 17:21

I just uploaded 2 files (CMakeLists.txt & main.cxx) which make it easy to reproduce the bug against a recent master in VTK 6 (efb1793c2215) and ITK v4.5.0.

Running the un-patched master I get the message "FAILED assert *pixelPtr == SecondFillValue" , with the patch applied there is no message.

(0032027)
Dave DeMarle (administrator)
2013-12-19 07:24

patch applied.
thanks John!

 Issue History
Date Modified Username Field Change
2011-10-28 12:02 John Stark New Issue
2011-10-28 12:07 John Stark Note Added: 0027620
2012-11-06 08:23 John Stark File Added: ComputePipelineMTime.patch
2012-11-06 08:31 John Stark Note Added: 0029655
2012-11-06 08:38 John Stark Note Edited: 0029655
2012-11-06 08:38 John Stark File Added: sampleCrash.cxx
2012-12-12 07:52 John Stark Note Edited: 0029655
2012-12-12 07:55 John Stark Note Added: 0029870
2013-12-16 12:57 Dave DeMarle Assigned To => Bill Lorensen
2013-12-16 12:57 Dave DeMarle Status backlog => tabled
2013-12-16 12:58 Dave DeMarle Status tabled => backlog
2013-12-16 16:41 John Stark Note Edited: 0029870
2013-12-18 17:13 John Stark File Added: CMakeLists.txt
2013-12-18 17:13 John Stark File Added: main.cxx
2013-12-18 17:16 John Stark Note Added: 0032026
2013-12-18 17:21 John Stark Note Edited: 0032026
2013-12-19 07:24 Dave DeMarle Note Added: 0032027
2013-12-19 07:24 Dave DeMarle Status backlog => closed
2013-12-19 07:24 Dave DeMarle Resolution open => fixed
2013-12-19 07:24 Dave DeMarle Fixed in Version => 6.1.0


Copyright © 2000 - 2018 MantisBT Team