[Paraview-developers] Bug found, then stuck...

Cornelis Bockemühl cornelis.bockemuehl at gmail.com
Sun Oct 22 09:25:27 EDT 2017


Dear VTK developers,

First of all: Am I posting this maybe in the wrong mailing list??

Anyway, here I go: In some project I am using the vtkIntersectPolyData
filter. This works as expected many times, but from time to time I had
crashes. I managed to track down these crashes to a little bug in that
class that can be fixed by a very little change:

in vtkIntersectionPolyDataFilter.cxx(14'1), replace

  lineBool = new bool[pd->GetNumberOfCells()];


by


  lineBool = new bool[pd->GetNumberOfCells() + 1];


The crash happened when this array was released (heap corruption), and it
turned out that indeed under some circumstances the program writes one
element more than there is allocated memory. This happens if further down
in the GetSingleLoop function the "else if" branch is chosen inside the "if
-- else if -- else" construct: Here the pd cells array is extended by one
element.

With the above fix the crash goes indeed away, but then later on the
program throws an exception at vtkCellLinks.h(171):

  /**

   * Increment the count of the number of cells using the point.

   */

  void IncrementLinkCount(vtkIdType ptId) { this->Array[ptId].ncells++;};


call stack at this point:


1   vtkCellLinks::IncrementLinkCount
vtkCellLinks.h                          171  0x7feced1c91a
2   vtkCellLinks::BuildLinks
vtkCellLinks.cxx                        153  0x7feced1a89a
3   vtkPolyData::BuildLinks
vtkPolyData.cxx                         1096 0x7fecf022ff1
4   debug_vtkIntersectionPolyDataFilter::Impl::SplitCell
debug_vtkIntersectionPolyDataFilter.cxx 1007 0x7fec0c994e2
5   debug_vtkIntersectionPolyDataFilter::Impl::SplitMesh
debug_vtkIntersectionPolyDataFilter.cxx 630  0x7fec0c971cc
6   debug_vtkIntersectionPolyDataFilter::RequestData
debug_vtkIntersectionPolyDataFilter.cxx 2590 0x7fec0c94531
7   vtkPolyDataAlgorithm::ProcessRequest
vtkPolyDataAlgorithm.cxx                90   0x7fecf556116
8   vtkExecutive::CallAlgorithm
vtkExecutive.cxx                        775  0x7fecf50eb55
9   vtkDemandDrivenPipeline::ExecuteData
vtkDemandDrivenPipeline.cxx             491  0x7fecf50477c
10  vtkCompositeDataPipeline::ExecuteData
vtkCompositeDataPipeline.cxx            179  0x7fecf4fc560
11  vtkDemandDrivenPipeline::ProcessRequest
vtkDemandDrivenPipeline.cxx             274  0x7fecf503343
12  vtkStreamingDemandDrivenPipeline::ProcessRequest
vtkStreamingDemandDrivenPipeline.cxx    328  0x7fecf588f67
13  vtkDemandDrivenPipeline::UpdateData
vtkDemandDrivenPipeline.cxx             444  0x7fecf5042da
14  vtkStreamingDemandDrivenPipeline::Update
vtkStreamingDemandDrivenPipeline.cxx    403  0x7fecf589412
15  vtkStreamingDemandDrivenPipeline::Update
vtkStreamingDemandDrivenPipeline.cxx    366  0x7fecf58916c
16  vtkAlgorithm::Update
vtkAlgorithm.cxx                        1457 0x7fecf4e5348
17  vtkAlgorithm::Update
vtkAlgorithm.cxx                        1451 0x7fecf4e53b3
18  vtkCemTestIntersectionFilter::RequestData
vtkCemTestIntersectionFilter.cxx        252  0x7fec0c84e7f
19  vtkPolyDataAlgorithm::ProcessRequest
vtkPolyDataAlgorithm.cxx                90   0x7fecf556116
20  vtkExecutive::CallAlgorithm
vtkExecutive.cxx                        775  0x7fecf50eb55
... <More>


Note that debug_vtkIntersectionPolyDataFilter is a copy of
vtkIntersectionPolyDataFilter that I prepared to do my testing and add some
debug output stuff etc. without interfering with the original code.

I guess that the exception is again a consequence of adding that "cell" in
GetSingleLoop, but now I must admit that my understanding of the internals
of these geometric object classes is not yet good enough to see what went
wrong! The code that needs checking by an expert seems to me the following

vtkIntersectionPolyDataFilter.cxx(1536..1555):

    //There is one line attached to point. This means the intersection has

    //an open intersection loop (i.e. the surfaces are open and one does not

    //completeley intersect the other.

    //Make an artificial triangle loop in this case

    else if (pointCells->GetNumberOfIds() < 2)

    {

      vtkSmartPointer<vtkPolyData> currentpd =

        vtkSmartPointer<vtkPolyData>::New();

      vtkSmartPointer<vtkCellArray> currentcells =

        vtkSmartPointer<vtkCellArray>::New();

      currentcells = pd->GetLines();

      currentcells->InsertNextCell(2);

      currentcells->InsertCellPoint(nextPt);

      currentcells->InsertCellPoint(startPt);

      nextCell = currentcells->GetNumberOfCells()-1;

      currentpd->SetLines(currentcells);

      currentpd->SetPoints(pd->GetPoints());

      pd->DeepCopy(currentpd);

      pd->BuildLinks();

    }


My guess is that this code is not yet good enough to properly insert the
"cell" into pd...

Any idea?

Regards,
Cornelis

-- 
Cornelis Bockemühl
Basel, Schweiz
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://public.kitware.com/pipermail/paraview-developers/attachments/20171022/6cb12e5a/attachment.html>


More information about the Paraview-developers mailing list