[vtkusers] vtkIntersectionPolyDataFilter::TriangleTriangleIntersection corrupts the stack?

Cory Quammen cquammen at cs.unc.edu
Mon Oct 8 16:39:36 EDT 2012


James,

This should be fixed in the master branch of VTK. Could you give it a try?

Thanks,
Cory

On Fri, Oct 5, 2012 at 3:48 PM, Cory Quammen <cquammen at cs.unc.edu> wrote:
> James,
>
> I've posted a patch on the gerrit review system here:
>
> http://review.source.kitware.com/#/c/7763/
>
> If you wouldn't mind, could you register on gerrit and give the patch a review?
>
> Additional instructions for VTK development that may be useful are posted here:
>
> http://www.vtk.org/Wiki/VTK/Git
>
> Thanks,
> Cory
>
> On Fri, Oct 5, 2012 at 11:41 AM, James Johnston <JamesJ at motionview3d.com> wrote:
>>> -----Original Message-----
>>> From: Cory Quammen [mailto:cquammen at cs.unc.edu]
>>> Sent: Thursday, October 04, 2012 20:19
>>> To: James Johnston
>>> Cc: vtkusers at vtk.org
>>> Subject: Re: [vtkusers]
>>> vtkIntersectionPolyDataFilter::TriangleTriangleIntersection corrupts the
>> stack?
>>>
>>> Hi James,
>>>
>>> I'm the author of this code. If you are interested in understanding the
>> algorithm
>>> more, please see Section 11.5.4 on pages 539-542 in
>>>
>>> Schnieder and Eberly, Geometric Tools for Computer Graphics, Morgan
>>> Kaufman
>>>
>>> I plotted your triangles in ParaView and they do not appear to intersect,
>> but
>>> rather share a vertex. Perhaps the line of intersection between the two
>> planes
>>> in which the triangles lie lines up with one of the edges. The three
>> edge/plane
>>> intersection tests could then all pass, causing overflow in the t1 and t2
>>> variables.
>>
>> Thanks for the response.  Apparently so.  They do indeed share a vertex and
>> I'm sure you noticed that the triangles are also nearly coplanar (or
>> coplanar I guess, depending on your precision).
>>
>>>
>>> I think the right thing to do would be to increase the size of t1 and
>>> t2 to 3 to avoid the stack overflow. Further down in the code, index1 and
>> index2
>>> are checked to see if they are equal to 2. If not, the triangle/triangle
>>> intersection would report that the triangles don't intersec, which seems
>>> reasonable to me.
>>>
>>> A more robust function would return a poly data that would contain a line
>>> segment for triangles that overlap nicely, a single vertex if the
>> triangles share a
>>> vertex, or several triangles if the triangles overlap and are co-planar.
>>> Unfortunately, I don't have the time to write that.
>>
>> Indeed, that would be most ideal.  Unfortunately, I don't have time to write
>> that for now at this point, either.  The main thing I'm focusing on finding
>> is the common case of triangles which intersect and don't share
>> edges/vertices and are not coplanar.  But it seems my data would make the
>> algorithm exhibit undefined behavior (stack corruption).
>>
>> As you indicated in your paper (which I skimmed through), there is a good
>> bit more future work to be done in this area for handling special cases like
>> coplanar triangles.
>>
>> For now, I would want the function to not crash and return *relatively* sane
>> results.  I say relatively since you already described the ideal.  For the
>> time being, I would think a patch to at least fix the crash would be
>> sufficient.
>>
>>>
>>> Thoughts?
>>>
>>> Cory
>>>
>>>
>>>
>>> On Tue, Oct 2, 2012 at 12:24 PM, James Johnston
>>> <JamesJ at motionview3d.com> wrote:
>>> > Hi,
>>> >
>>> > I'm trying to find triangle-triangle intersection and I think
>>> > vtkIntersectionPolyDataFilter::TriangleTriangleIntersection will work
>>> > for me.  Unfortunately, I *think* this function contains a bug.  I
>>> > just don't know the algorithms in that function well enough to say
>>> > what the exact problem is or what the fix should be - I'm just trying to
>> use the
>>> function.
>>> > I assume there aren't any bugs in my code - but if there are,
>>> > hopefully someone can point them out!
>>> >
>>> > Trivial code that causes crash:
>>> >
>>> >         int coplanar;
>>> >         double isectpt1[3], isectpt2[3];
>>> >         double thisCellTri[9] = {-30.121610641479492,
>>> > 29.366317749023438, -27.164772033691406, -29.923696517944336,
>>> > 29.366317749023438, -27.334911346435547, -30.055425643920898,
>>> > 28.556877136230469, -27.227336883544922};
>>> >         double otherCellTri[9] = {-29.923696517944336,
>>> > 29.366317749023438, -27.334911346435547, -29.858596801757812,
>>> > 29.366317749023438, -27.390081405639648, -29.813289642333984,
>>> > 27.682807922363281, -27.436887741088867};
>>> >
>>> > vtkIntersectionPolyDataFilter::TriangleTriangleIntersection(&thisCellT
>>> > ri[0], &thisCellTri[3], &thisCellTri[6], &otherCellTri[0],
>>> > &otherCellTri[3], &otherCellTri[6],
>>> >                             coplanar, isectpt1, isectpt2);
>>> >
>>> > This causes TriangleTriangleIntersection to corrupt the stack.  The
>>> > offending section of code inside TriangleTriangleIntersection is as
>> follows:
>>> >
>>> >   int index1 = 0, index2 = 0;
>>> >   double t1[2], t2[2];
>>> >   for (int i = 0; i < 3; i++)
>>> >     {
>>> >     double t, x[3];
>>> >     int id1 = i, id2 = (i+1) % 3;
>>> >
>>> >     // Find t coordinate on line of intersection between two planes.
>>> >     if (vtkPlane::IntersectWithLine( pts1[id1], pts1[id2], n2, p2, t, x
>> ))
>>> >       {
>>> >       t1[index1++] = vtkMath::Dot(x, v) - vtkMath::Dot(p, v);
>>> >       }
>>> >
>>> >     if (vtkPlane::IntersectWithLine( pts2[id1], pts2[id2], n1, p1, t, x
>> ))
>>> >       {
>>> >       t2[index2++] = vtkMath::Dot(x, v) - vtkMath::Dot(p, v);
>>> >       }
>>> >     }
>>> >
>>> > I've verified that the first " vtkPlane::IntersectWithLine( pts1[id1],
>>> > pts1[id2], n2, p2, t, x )" condition is evaluating to true all three
>> times.
>>> > This causes the loop to attempt to write to t1[2], which is
>>> > out-of-bounds and therefore corrupts the stack (detected by MSVC debug
>>> > build).  I don't want to just change the declaration of t1 and t2 to
>>> > hold 3 elements instead of 2, because the rest of the code doesn't
>>> > seem to reference a third element.  It would fix the stack corruption
>>> > but probably just mask the real issue.
>>> >
>>> > This is tested on VTK 5.10.0.  However, I assume the latest version on
>>> > master branch contains the same issue, because Git history shows that
>>> > nobody has made any meaningful changes to
>>> > vtkIntersectionPolyDataFilter.cxx.  (Just a few changes for VTK6 and
>>> > to fix an implicit double-to-int compiler warning).
>>> >
>>> > Ideas?
>>> >
>>> > Best regards,
>>> >
>>> > James Johnston
>>
>> _______________________________________________
>> Powered by www.kitware.com
>>
>> Visit other Kitware open-source projects at http://www.kitware.com/opensource/opensource.html
>>
>> Please keep messages on-topic and check the VTK FAQ at: http://www.vtk.org/Wiki/VTK_FAQ
>>
>> Follow this link to subscribe/unsubscribe:
>> http://www.vtk.org/mailman/listinfo/vtkusers
>
>
>
> --
> Cory Quammen
> Research Associate
> Department of Computer Science
> The University of North Carolina at Chapel Hill



-- 
Cory Quammen
Research Associate
Department of Computer Science
The University of North Carolina at Chapel Hill



More information about the vtkusers mailing list