Skip to content
Snippets Groups Projects
Commit 761b5700 authored by David Gobbi's avatar David Gobbi
Browse files

Fix streaming in vtkImageMarchingCubes

This filter wasn't streaming because the input's whole extent was
updated before vtkImageMarchingCubes::RequestData() was called.
In order for streaming to work properly, the input data must not
be updated before RequestData() updates it.  To state this in
another way, setting a memory limit for the input has no effect
if the whole input is filled before vtkImageMarchingCubes even
executes.

To properly test this filter, the test was modified to use a
streaming reader (vtkImageReader2 instead of vtkVolume16Reader),
and the InputMemoryLimit was reduced.
parent 8c78f17d
Branches
Tags
No related merge requests found
......@@ -70,6 +70,21 @@ int vtkImageMarchingCubesGetTypeSize(T*)
return sizeof(T);
}
//----------------------------------------------------------------------------
int vtkImageMarchingCubes::RequestUpdateExtent(
vtkInformation* vtkNotUsed(request),
vtkInformationVector** inputVector,
vtkInformationVector* vtkNotUsed(outputVector))
{
// start with an empty UPDATE_EXTENT to ensure proper streaming
// (the UPDATE_EXTENT will be set properly in RequestData()).
int extent[6] = { 0, -1, 0, -1, 0, -1 };
vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_EXTENT(),
extent, 6);
return 1;
}
//----------------------------------------------------------------------------
int vtkImageMarchingCubes::RequestData(
vtkInformation *vtkNotUsed(request),
......@@ -123,6 +138,7 @@ int vtkImageMarchingCubes::RequestData(
return 1;
}
// Get the WHOLE_EXTENT and divide it into chunks
int extent[6];
inInfo->Get(vtkStreamingDemandDrivenPipeline::WHOLE_EXTENT(), extent);
// multiply by the area of each slice
......
......@@ -146,7 +146,10 @@ protected:
int LocatorMinX;
int LocatorMinY;
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override;
int RequestData(vtkInformation *, vtkInformationVector **,
vtkInformationVector *) override;
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **,
vtkInformationVector *) override;
int FillInputPortInformation(int port, vtkInformation *info) override;
void March(vtkImageData *inData, int chunkMin, int chunkMax,
......
......@@ -22,19 +22,17 @@ iren.SetRenderWindow(renWin)
# create pipeline
#
v16 = vtk.vtkVolume16Reader()
v16.SetDataDimensions(64, 64)
v16.GetOutput().SetOrigin(0.0, 0.0, 0.0)
v16.SetDataByteOrderToLittleEndian()
v16.SetFilePrefix(VTK_DATA_ROOT + "/Data/headsq/quarter")
v16.SetImageRange(1, 93)
v16.SetDataSpacing(3.2, 3.2, 1.5)
v16.Update()
reader = vtk.vtkImageReader2()
reader.SetDataScalarTypeToUnsignedShort()
reader.SetFilePrefix(VTK_DATA_ROOT + "/Data/headsq/quarter")
reader.SetDataExtent(0, 63, 0, 63, 1, 93)
reader.SetDataSpacing(3.2, 3.2, 1.5)
reader.SetDataOrigin(0.0, 0.0, 0.0)
iso = vtk.vtkImageMarchingCubes()
iso.SetInputConnection(v16.GetOutputPort())
iso.SetInputConnection(reader.GetOutputPort())
iso.SetValue(0, 1150)
iso.SetInputMemoryLimit(1000)
iso.SetInputMemoryLimit(100)
isoMapper = vtk.vtkPolyDataMapper()
isoMapper.SetInputConnection(iso.GetOutputPort())
......@@ -45,7 +43,7 @@ isoActor.SetMapper(isoMapper)
isoActor.GetProperty().SetColor(GetRGBColor('antique_white'))
outline = vtk.vtkOutlineFilter()
outline.SetInputConnection(v16.GetOutputPort())
outline.SetInputConnection(reader.GetOutputPort())
outlineMapper = vtk.vtkPolyDataMapper()
outlineMapper.SetInputConnection(outline.GetOutputPort())
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment