Skip to content
Snippets Groups Projects
Commit 3df199ea authored by Andrew Maclean's avatar Andrew Maclean Committed by Code Review
Browse files

Merge topic '0014457-ExtractPolyData-PassPoints' into master

0567a19c Extended filter to support culling unused points.
parents dd0222fb 0567a19c
No related branches found
No related tags found
No related merge requests found
a97bdd5c11cf4b60d4ecb824e7771a4a
1a0557264d4a8511a547dd5744bcdcb3
......@@ -12,6 +12,7 @@ sphere = vtk.vtkSphereSource()
sphere.SetThetaResolution(8)
sphere.SetPhiResolution(16)
sphere.SetRadius(1.5)
# Extraction stuff
t = vtk.vtkTransform()
t.RotateX(90)
......@@ -22,11 +23,48 @@ extract = vtk.vtkExtractPolyDataGeometry()
extract.SetInputConnection(sphere.GetOutputPort())
extract.SetImplicitFunction(cylfunc)
extract.ExtractBoundaryCellsOn()
extract.PassPointsOn()
sphereMapper = vtk.vtkPolyDataMapper()
sphereMapper.SetInputConnection(extract.GetOutputPort())
sphereMapper.GlobalImmediateModeRenderingOn()
sphereActor = vtk.vtkActor()
sphereActor.SetMapper(sphereMapper)
# Extraction stuff - now cull points
extract2 = vtk.vtkExtractPolyDataGeometry()
extract2.SetInputConnection(sphere.GetOutputPort())
extract2.SetImplicitFunction(cylfunc)
extract2.ExtractBoundaryCellsOn()
extract2.PassPointsOff()
sphereMapper2 = vtk.vtkPolyDataMapper()
sphereMapper2.SetInputConnection(extract2.GetOutputPort())
sphereActor2 = vtk.vtkActor ()
sphereActor2.SetMapper(sphereMapper2)
sphereActor2.AddPosition(2.5, 0, 0)
# Put some glyphs on the points
glyphSphere = vtk.vtkSphereSource()
glyphSphere.SetRadius(0.05)
glyph = vtk.vtkGlyph3D()
glyph.SetInputConnection(extract.GetOutputPort())
glyph.SetSourceConnection(glyphSphere.GetOutputPort())
glyph.SetScaleModeToDataScalingOff()
glyphMapper = vtk.vtkPolyDataMapper()
glyphMapper.SetInputConnection(glyph.GetOutputPort())
glyphActor = vtk.vtkActor()
glyphActor.SetMapper(glyphMapper)
glyph2 = vtk.vtkGlyph3D()
glyph2.SetInputConnection(extract2.GetOutputPort())
glyph2.SetSourceConnection(glyphSphere.GetOutputPort())
glyph2.SetScaleModeToDataScalingOff()
glyphMapper2 = vtk.vtkPolyDataMapper()
glyphMapper2.SetInputConnection(glyph2.GetOutputPort())
glyphActor2 = vtk.vtkActor()
glyphActor2.SetMapper(glyphMapper2)
glyphActor2.AddPosition(2.5, 0, 0)
# Create the RenderWindow, Renderer and both Actors
#
ren1 = vtk.vtkRenderer()
......@@ -38,6 +76,9 @@ iren.SetRenderWindow(renWin)
# Add the actors to the renderer, set the background and size
#
ren1.AddActor(sphereActor)
ren1.AddActor(glyphActor)
ren1.AddActor(sphereActor2)
ren1.AddActor(glyphActor2)
ren1.ResetCamera()
ren1.GetActiveCamera().Azimuth(30)
ren1.SetBackground(0.1,0.2,0.4)
......
......@@ -11,7 +11,7 @@ vtkSphereSource sphere
sphere SetPhiResolution 16
sphere SetRadius 1.5
# Extraction stuff
# Extraction stuff - first just pass points
vtkTransform t
t RotateX 90
vtkCylinder cylfunc
......@@ -21,14 +21,47 @@ vtkExtractPolyDataGeometry extract
extract SetInputConnection [sphere GetOutputPort]
extract SetImplicitFunction cylfunc
extract ExtractBoundaryCellsOn
extract PassPointsOn
vtkPolyDataMapper sphereMapper
sphereMapper SetInputConnection [extract GetOutputPort]
sphereMapper GlobalImmediateModeRenderingOn
vtkActor sphereActor
sphereActor SetMapper sphereMapper
# Extraction stuff - now cull points
vtkExtractPolyDataGeometry extract2
extract2 SetInputConnection [sphere GetOutputPort]
extract2 SetImplicitFunction cylfunc
extract2 ExtractBoundaryCellsOn
extract2 PassPointsOff
vtkPolyDataMapper sphereMapper2
sphereMapper2 SetInputConnection [extract2 GetOutputPort]
vtkActor sphereActor2
sphereActor2 SetMapper sphereMapper2
sphereActor2 AddPosition 2.5 0 0
# Put some glyphs on the points
vtkSphereSource glyphSphere
glyphSphere SetRadius 0.05
vtkGlyph3D glyph
glyph SetInputConnection [extract GetOutputPort]
glyph SetSourceConnection [glyphSphere GetOutputPort]
glyph SetScaleModeToDataScalingOff
vtkPolyDataMapper glyphMapper
glyphMapper SetInputConnection [glyph GetOutputPort]
vtkActor glyphActor
glyphActor SetMapper glyphMapper
vtkGlyph3D glyph2
glyph2 SetInputConnection [extract2 GetOutputPort]
glyph2 SetSourceConnection [glyphSphere GetOutputPort]
glyph2 SetScaleModeToDataScalingOff
vtkPolyDataMapper glyphMapper2
glyphMapper2 SetInputConnection [glyph2 GetOutputPort]
vtkActor glyphActor2
glyphActor2 SetMapper glyphMapper2
glyphActor2 AddPosition 2.5 0 0
# Create the RenderWindow, Renderer and both Actors
#
vtkRenderer ren1
......@@ -41,6 +74,9 @@ vtkRenderWindowInteractor iren
# Add the actors to the renderer, set the background and size
#
ren1 AddActor sphereActor
ren1 AddActor glyphActor
ren1 AddActor sphereActor2
ren1 AddActor glyphActor2
ren1 ResetCamera
[ren1 GetActiveCamera] Azimuth 30
......
......@@ -28,6 +28,7 @@ vtkStandardNewMacro(vtkExtractPolyDataGeometry);
vtkCxxSetObjectMacro(vtkExtractPolyDataGeometry,
ImplicitFunction,vtkImplicitFunction);
//----------------------------------------------------------------------------
// Construct object with ExtractInside turned on.
vtkExtractPolyDataGeometry::vtkExtractPolyDataGeometry(vtkImplicitFunction *f)
{
......@@ -39,13 +40,16 @@ vtkExtractPolyDataGeometry::vtkExtractPolyDataGeometry(vtkImplicitFunction *f)
this->ExtractInside = 1;
this->ExtractBoundaryCells = 0;
this->PassPoints = 0;
}
//----------------------------------------------------------------------------
vtkExtractPolyDataGeometry::~vtkExtractPolyDataGeometry()
{
this->SetImplicitFunction(NULL);
}
//----------------------------------------------------------------------------
// Overload standard modified time function. If implicit function is modified,
// then this object is modified as well.
unsigned long vtkExtractPolyDataGeometry::GetMTime()
......@@ -62,6 +66,7 @@ unsigned long vtkExtractPolyDataGeometry::GetMTime()
return mTime;
}
//----------------------------------------------------------------------------
int vtkExtractPolyDataGeometry::RequestData(
vtkInformation *vtkNotUsed(request),
vtkInformationVector **inputVector,
......@@ -82,10 +87,11 @@ int vtkExtractPolyDataGeometry::RequestData(
vtkPointData *outputPD = output->GetPointData();
vtkCellData *outputCD = output->GetCellData();
vtkPoints *inPts=input->GetPoints();
vtkIdType numPts, i, cellId = 0, newId;
vtkIdType numPts, i, cellId = 0, newId, ptId, *pointMap=NULL;
float multiplier;
vtkCellArray *inVerts=NULL, *inLines=NULL, *inPolys=NULL, *inStrips=NULL;
vtkCellArray *newVerts=NULL, *newLines=NULL, *newPolys=NULL, *newStrips=NULL;
vtkPoints *newPts=NULL;
vtkDebugMacro(<< "Extracting poly data geometry");
......@@ -111,14 +117,35 @@ int vtkExtractPolyDataGeometry::RequestData(
vtkFloatArray *newScalars = vtkFloatArray::New();
newScalars->SetNumberOfValues(numPts);
for (int ptId=0; ptId < numPts; ptId++ )
for (ptId=0; ptId < numPts; ptId++ )
{
newScalars->SetValue(ptId, this->ImplicitFunction->
FunctionValue(inPts->GetPoint(ptId))*multiplier);
}
output->SetPoints(inPts);
outputPD->PassData(pd);
// Do different things with the points depending on user directive
if ( this->PassPoints )
{
output->SetPoints(inPts);
outputPD->PassData(pd);
}
else
{
newPts = vtkPoints::New();
newPts->Allocate(numPts/4,numPts);
pointMap = new vtkIdType[numPts]; // maps old point ids into new
for (ptId=0; ptId < numPts; ptId++)
{
if ( newScalars->GetValue(ptId) <= 0.0 )
{
newId = this->InsertPointInMap(ptId, inPts, newPts, pointMap);
}
else
{
pointMap[ptId] = -1;
}
}
}
outputCD->CopyAllocate(cd);
// Now loop over all cells to see whether they are inside the implicit
......@@ -168,7 +195,26 @@ int vtkExtractPolyDataGeometry::RequestData(
}
if ( (numIn == npts) || (this->ExtractBoundaryCells && numIn > 0) )
{
newId = newVerts->InsertNextCell(npts,pts);
if ( this->PassPoints )
{
newId = newVerts->InsertNextCell(npts,pts);
}
else
{
newId = newVerts->InsertNextCell(npts);
for ( i=0; i < npts; i++)
{
if ( pointMap[pts[i]] < 0 )
{
ptId = this->InsertPointInMap(pts[i], inPts, newPts, pointMap);
}
else
{
ptId = pointMap[pts[i]];
}
newVerts->InsertCellPoint(ptId);
}
}
outputCD->CopyData(cd, cellId, newId);
}
cellId++;
......@@ -190,7 +236,26 @@ int vtkExtractPolyDataGeometry::RequestData(
}
if ( (numIn == npts) || (this->ExtractBoundaryCells && numIn > 0) )
{
newId = newLines->InsertNextCell(npts,pts);
if ( this->PassPoints )
{
newId = newLines->InsertNextCell(npts,pts);
}
else
{
newId = newLines->InsertNextCell(npts);
for ( i=0; i < npts; i++)
{
if ( pointMap[pts[i]] < 0 )
{
ptId = this->InsertPointInMap(pts[i], inPts, newPts, pointMap);
}
else
{
ptId = pointMap[pts[i]];
}
newLines->InsertCellPoint(ptId);
}
}
outputCD->CopyData(cd, cellId, newId);
}
cellId++;
......@@ -212,7 +277,26 @@ int vtkExtractPolyDataGeometry::RequestData(
}
if ( (numIn == npts) || (this->ExtractBoundaryCells && numIn > 0) )
{
newId = newPolys->InsertNextCell(npts,pts);
if ( this->PassPoints )
{
newId = newPolys->InsertNextCell(npts,pts);
}
else
{
newId = newPolys->InsertNextCell(npts);
for ( i=0; i < npts; i++)
{
if ( pointMap[pts[i]] < 0 )
{
ptId = this->InsertPointInMap(pts[i], inPts, newPts, pointMap);
}
else
{
ptId = pointMap[pts[i]];
}
newPolys->InsertCellPoint(ptId);
}
}
outputCD->CopyData(cd, cellId, newId);
}
cellId++;
......@@ -234,7 +318,26 @@ int vtkExtractPolyDataGeometry::RequestData(
}
if ( (numIn == npts) || (this->ExtractBoundaryCells && numIn > 0) )
{
newId = newStrips->InsertNextCell(npts,pts);
if ( this->PassPoints )
{
newId = newStrips->InsertNextCell(npts,pts);
}
else
{
newId = newStrips->InsertNextCell(npts);
for ( i=0; i < npts; i++)
{
if ( pointMap[pts[i]] < 0 )
{
ptId = this->InsertPointInMap(pts[i], inPts, newPts, pointMap);
}
else
{
ptId = pointMap[pts[i]];
}
newStrips->InsertCellPoint(ptId);
}
}
outputCD->CopyData(cd, cellId, newId);
}
cellId++;
......@@ -245,6 +348,20 @@ int vtkExtractPolyDataGeometry::RequestData(
// Update ourselves and release memory
//
newScalars->Delete();
if ( ! this->PassPoints )
{
output->SetPoints(newPts);
newPts->Delete();
outputPD->CopyAllocate(pd);
for (i=0; i < numPts; i++)
{
if ( pointMap[i] >= 0 )
{
outputPD->CopyData(pd,i,pointMap[i]);
}
}
delete [] pointMap;
}
if ( newVerts )
{
......@@ -270,6 +387,7 @@ int vtkExtractPolyDataGeometry::RequestData(
return 1;
}
//----------------------------------------------------------------------------
void vtkExtractPolyDataGeometry::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os,indent);
......@@ -287,4 +405,6 @@ void vtkExtractPolyDataGeometry::PrintSelf(ostream& os, vtkIndent indent)
<< (this->ExtractInside ? "On\n" : "Off\n");
os << indent << "Extract Boundary Cells: "
<< (this->ExtractBoundaryCells ? "On\n" : "Off\n");
os << indent << "Pass Points: "
<< (this->PassPoints ? "On\n" : "Off\n");
}
......@@ -26,6 +26,11 @@
// region.) An option exists to extract cells that are neither inside nor
// outside (i.e., boundary).
//
// Note that this filter also has the option to directly pass all points or cull
// the points that do not satisfy the implicit function test. Passing all points
// is a tad faster, but then points remain that do not pass the test and may mess
// up subsequent glyphing operations and so on. By default points are culled.
//
// A more general version of this filter is available for arbitrary
// vtkDataSet input (see vtkExtractGeometry).
......@@ -74,6 +79,13 @@ public:
vtkGetMacro(ExtractBoundaryCells,int);
vtkBooleanMacro(ExtractBoundaryCells,int);
// Description:
// Boolean controls whether points are culled or simply passed through
// to the output.
vtkSetMacro(PassPoints,int);
vtkGetMacro(PassPoints,int);
vtkBooleanMacro(PassPoints,int);
protected:
vtkExtractPolyDataGeometry(vtkImplicitFunction *f=NULL);
~vtkExtractPolyDataGeometry();
......@@ -84,9 +96,25 @@ protected:
vtkImplicitFunction *ImplicitFunction;
int ExtractInside;
int ExtractBoundaryCells;
int PassPoints;
vtkIdType InsertPointInMap(vtkIdType i, vtkPoints *inPts, vtkPoints *newPts, vtkIdType *pointMap);
private:
vtkExtractPolyDataGeometry(const vtkExtractPolyDataGeometry&); // Not implemented.
void operator=(const vtkExtractPolyDataGeometry&); // Not implemented.
};
// Description:
// When not passing points, have to use a point map to keep track of things.
inline vtkIdType vtkExtractPolyDataGeometry::InsertPointInMap(vtkIdType i, vtkPoints *inPts,
vtkPoints *newPts, vtkIdType *pointMap)
{
double x[3];
inPts->GetPoint(i, x);
pointMap[i] = newPts->InsertNextPoint(x);
return pointMap[i];
}
#endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment