[vtk-developers] bug in Marching Cubes
Eduardo Suarez-Santana
esuarez at itccanarias.org
Wed Mar 30 08:14:10 EDT 2011
El 30/03/11 12:46, David Doria escribió:
> On Wed, Mar 30, 2011 at 6:48 AM, Eduardo Suarez-Santana
> <esuarez at itccanarias.org <mailto:esuarez at itccanarias.org>> wrote:
>
> In vtk-5.6.1/Graphics/MarchingCubes.cxx:
>
> 499 vtkDoubleArray *image=vtkDoubleArray::New();
> 500
> image->SetNumberOfComponents(inScalars->GetNumberOfComponents());
> 501
> image->SetNumberOfTuples(image->GetNumberOfComponents()*dataSize);
> 502 inScalars->GetTuples(0,dataSize,image);
> 503
> 504 double *scalars = image->GetPointer(0);
>
> I think it should be:
>
> 499 vtkDoubleArray *image=vtkDoubleArray::New();
> 500
> image->SetNumberOfComponents(inScalars->GetNumberOfComponents());
> 501 image->SetNumberOfTuples(dataSize); //changed
> 502 inScalars->GetTuples(0,dataSize-1,image); // changed
> 503
> 504 double *scalars = image->GetPointer(0);
>
>
> Why do you think that? Do you have a test that fails the old way and
> passes the new way?
>
> David
IMHO it is a "semantic" bug.
It makes no sense because:
1) It reserves too much memory (501).
dataSize*numComps tuples, that is, dataSize*numComps*numComps elements.
However,
dataSize tuples, that is, dataSize*numComps elements is enough.
2) it is getting tuples beyond inScalars range (502).
inScalars are the Scalars of the vtkImageData input, of size dataSize,
so dataSize tuples are expected. Accessing vtkIdType number dataSize is
out of range.
Moreover:
3) It does not fail the new way, at least in my cases.
And anyway,
4) (504) is providing weird input to marching cubes.
GetPointer(0), for a multiple component image, provides a pointer to a
mix of components and tuples. Getting the first dataSize elements of the
array (coded in variable dims, call in line 505) is going to mix
components and tuples making a nosense input:
505
vtkMarchingCubesComputeGradient(this,scalars,dims,origin,spacing,this->Locator,
506 newScalars,newGradients,
507 newNormals,newPolys,values,numContours);
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://public.kitware.com/pipermail/vtk-developers/attachments/20110330/2b13c273/attachment.html>
More information about the vtk-developers
mailing list