[vtkusers] Contours and vtkFeatureEdges

johnny carroll norris jnorris at csar.uiuc.edu
Fri May 25 18:12:33 EDT 2001


Hi,

One of the functions of my vis package is to generate isosurfaces with
"highlights," which are basically just the isosurfaces plus tubes wrapped
around the isosurface edges.  My problem is that for cell types other than
vtkHexahedron, the Contour() method generates extraneous edges.  This creates
a "shattered glass" effect when generating isosurfaces with highlights.  The
first attached image illustrates this problem.  The vtkPyramid and vtkWedge
cells have the same problem, but I don't have an example for these.

Since the problem doesn't exist for vtkHexahedron, I patched
vtkTetra::Contour() to use the same technique at vtkHexahedron::Contour().
This fixes the problem with the extra edges, but there's now a coloring
problem.  The second attached image illustrates.

I'm also attaching the diff of the original vtkTetra.cxx and my patched
version, as well as the example program I used to generate these images.

I'm using nightly release of VTK3.2 from May 9th, Solaris 7, and CC from
SC 5.0.  Any suggestions would be greatly appreciated.

Thanks,
John

P.S. My patch also fixes a problem I reported months ago, where the triangles
generated by vtkTetra::Contour() face the opposite way of triangles generated
by the other cells' Contour() methods.  This was a problem when using
vtkCutter on a vtkUnstructuredGrid with a mixture of cell types, and setting
BackfaceCullingOn() in the actor's properties.
-- 
John Norris
Research Programmer
Center for Simulation of Advanced Rockets
http://www.uiuc.edu/ph/www/jnorris
-------------- next part --------------
A non-text attachment was scrubbed...
Name: bad.gif
Type: image/gif
Size: 11951 bytes
Desc: not available
URL: <http://www.vtk.org/pipermail/vtkusers/attachments/20010525/0cc042cf/attachment.gif>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: better.gif
Type: image/gif
Size: 15640 bytes
Desc: not available
URL: <http://www.vtk.org/pipermail/vtkusers/attachments/20010525/0cc042cf/attachment-0001.gif>
-------------- next part --------------
*** vtkTetra.cxx.orig	Fri May 25 16:38:55 2001
--- vtkTetra.cxx	Fri May 25 16:43:22 2001
***************
*** 297,305 ****
    static int CASE_MASK[4] = {1,2,4,8};
    TRIANGLE_CASES *triCase;
    EDGE_LIST  *edge;
!   int i, j, index, *vert, newCellId;
    int pts[3];
!   float t, *x1, *x2, x[3];
  
    // Build the case table
    for ( i=0, index = 0; i < 4; i++)
--- 297,305 ----
    static int CASE_MASK[4] = {1,2,4,8};
    TRIANGLE_CASES *triCase;
    EDGE_LIST  *edge;
!   int i, j, k, index, *vert, e1, e2, newCellId;
    int pts[3];
!   float t, x1[3], x2[3], x[3], deltaScalar;
  
    // Build the case table
    for ( i=0, index = 0; i < 4; i++)
***************
*** 315,338 ****
  
    for ( ; edge[0] > -1; edge += 3 )
      {
!     for (i=0; i<3; i++) // insert triangle
        {
        vert = edges[edge[i]];
!       t = (value - cellScalars->GetScalar(vert[0])) /
!           (cellScalars->GetScalar(vert[1]) - cellScalars->GetScalar(vert[0]));
!       x1 = this->Points->GetPoint(vert[0]);
!       x2 = this->Points->GetPoint(vert[1]);
        for (j=0; j<3; j++)
  	{
  	x[j] = x1[j] + t * (x2[j] - x1[j]);
  	}
!       if ( locator->InsertUniquePoint(x, pts[i]) )
          {
          if ( outPd ) 
            {
!           int p1 = this->PointIds->GetId(vert[0]);
!           int p2 = this->PointIds->GetId(vert[1]);
!           outPd->InterpolateEdge(inPd,pts[i],p1,p2,t);
            }
          }
        }
--- 315,359 ----
  
    for ( ; edge[0] > -1; edge += 3 )
      {
!     for (i=0,k=2; i<3; i++,k--) // insert triangle
        {
        vert = edges[edge[i]];
!       // calculate a preferred interpolation direction
!       deltaScalar = cellScalars->GetScalar(vert[1]) - cellScalars->GetScalar(vert[0]);
!       if (deltaScalar > 0)
!         {
!         e1 = vert[0]; e2 = vert[1];
!         }
!       else
!         {
!         e1 = vert[1]; e2 = vert[0];
!         deltaScalar = -deltaScalar;
!         }
! 
!       // linear interpolation
!       if (deltaScalar == 0.0)
!         {
!         t = 0.0;
!         }
!       else
!         {
!         t = (value - cellScalars->GetScalar(e1)) / deltaScalar;
!         }
! 
!       this->Points->GetPoint(e1, x1);
!       this->Points->GetPoint(e2, x2);
! 
        for (j=0; j<3; j++)
  	{
  	x[j] = x1[j] + t * (x2[j] - x1[j]);
  	}
!       if ( locator->InsertUniquePoint(x, pts[k]) )
          {
          if ( outPd ) 
            {
!           int p1 = this->PointIds->GetId(vert[e1]);
!           int p2 = this->PointIds->GetId(vert[e2]);
!           outPd->InterpolateEdge(inPd,pts[k],p1,p2,t);
            }
          }
        }
-------------- next part --------------
#include <vtkSphereSource.h>
#include <vtkPolyData.h>
#include <vtkPoints.h>
#include <vtkUnstructuredGrid.h>
#include <vtkIdList.h>
#include <vtkSubdivideTetra.h>
#include <vtkContourGrid.h>
#include <vtkFeatureEdges.h>
#include <vtkTubeFilter.h>
#include <vtkAppendPolyData.h>
#include <vtkPolyDataMapper.h>
#include <vtkActor.h>
#include <vtkRenderer.h>
#include <vtkRenderWindow.h>
#include <vtkRenderWindowInteractor.h>

int main()
{
   vtkSphereSource *pSphere = vtkSphereSource::New();
   pSphere->SetRadius(50);

   vtkPolyData *pData = pSphere->GetOutput();
   pData->Update();

   vtkPoints *pPoints = vtkPoints::New();
   pPoints->DeepCopy(pData->GetPoints());
   pPoints->InsertNextPoint(0.0, 0.0, 0.0);

   vtkUnstructuredGrid *pGrid = vtkUnstructuredGrid::New();
   pGrid->Allocate(pData->GetNumberOfCells());
   pGrid->SetPoints(pPoints);
   pPoints->Delete();

   int i;
   vtkIdList *pIds = vtkIdList::New();
   for (i=0; i<pData->GetNumberOfCells(); ++i) {
      if (pData->GetCellType(i) != VTK_TRIANGLE) {
         vtkGenericWarningMacro(<< "Cell type is " << pData->GetCellType(i));
         pIds->Delete();
         pGrid->Delete();
         pSphere->Delete();
         return 1;
      }

      pData->GetCellPoints(i, pIds);
      pIds->InsertNextId(0);

      pGrid->InsertNextCell(VTK_TETRA, pIds);
   }
   pIds->Delete();

   vtkSubdivideTetra *pDiv1 = vtkSubdivideTetra::New();
   pDiv1->SetInput(pGrid);
   pGrid->Delete();

   vtkSubdivideTetra *pDiv2 = vtkSubdivideTetra::New();
   pDiv2->SetInput(pDiv1->GetOutput());
   pDiv1->Delete();

   pGrid = pDiv2->GetOutput();
   pGrid->Update();

   pPoints = pGrid->GetPoints();
   vtkScalars *pScalars = vtkScalars::New();
   pScalars->SetNumberOfScalars(pPoints->GetNumberOfPoints());
   for (i=0; i<pPoints->GetNumberOfPoints(); ++i)
      pScalars->SetScalar(i, pPoints->GetPoint(i)[0]);

   pGrid->GetPointData()->SetScalars(pScalars);

   vtkContourGrid *pContour = vtkContourGrid::New();
   pContour->GenerateValues(6, pScalars->GetRange());
   pContour->SetInput(pGrid);
   pDiv2->Delete();

   vtkFeatureEdges *pEdges = vtkFeatureEdges::New();
   pEdges->BoundaryEdgesOn();
   pEdges->FeatureEdgesOff();
   pEdges->NonManifoldEdgesOff();
   pEdges->ManifoldEdgesOff();
   pEdges->SetInput(pContour->GetOutput());

   vtkTubeFilter *pTube = vtkTubeFilter::New();
   pTube->SetNumberOfSides(8);
   pTube->SetRadius(1.f);
   pTube->SetVaryRadiusToVaryRadiusOff();
   pTube->SetInput(pEdges->GetOutput());
   pEdges->Delete();

   vtkAppendPolyData *pAppend = vtkAppendPolyData::New();
   pAppend->AddInput(pContour->GetOutput());
   pAppend->AddInput(pTube->GetOutput());
   pContour->Delete();
   pTube->Delete();

   vtkPolyDataMapper *pMapper = vtkPolyDataMapper::New();
   pMapper->SetScalarRange(pScalars->GetRange());
   pMapper->SetInput(pAppend->GetOutput());
   pScalars->Delete();
   pAppend->Delete();

   vtkActor *pActor = vtkActor::New();
   pActor->SetMapper(pMapper);
   pMapper->Delete();

   vtkRenderer *pRenderer = vtkRenderer::New();
   pRenderer->AddActor(pActor);
   pActor->Delete();

   vtkRenderWindow *pWindow = vtkRenderWindow::New();
   pWindow->AddRenderer(pRenderer);
   pRenderer->Delete();

   vtkRenderWindowInteractor *pIRen = vtkRenderWindowInteractor::New();
   pIRen->SetRenderWindow(pWindow);

   pWindow->Render();
   pWindow->Delete();

   pIRen->Start();

   pIRen->Delete();
   pSphere->Delete();

   return 0;
}


More information about the vtkusers mailing list