Skip to content

ClipUnstructuredGridWithPlane

vtk-examples/Python/UnstructuredGrid/ClipUnstructuredGridWithPlane

Description

The example uses vtkTableBasedClipDataSet to clip a vtkUnstructuredGrid. The resulting output and clipped output are presented in yellow and red respectively. To illustrate the clipped interfaces, the example uses a vtkTransform to rotate each output about their centers.

Note that unlike other clipping filters (except for vtkClipPolyData), vtkTableBasedClipDataSet retains the original cells if they are not clipped.

After exiting, the example reports the number of each cell type for each output:

The inside dataset contains a
  [vtkUnstructuredGrid](https://www.vtk.org/doc/nightly/html/classvtkUnstructuredGrid.html#details) that has 26116 cells
    Cell type [vtkTetra](https://www.vtk.org/doc/nightly/html/classvtkTetra.html#details) occurs 3751 times.
    Cell type [vtkHexahedron](https://www.vtk.org/doc/nightly/html/classvtkHexahedron.html#details) occurs 17361 times.
    Cell type [vtkWedge](https://www.vtk.org/doc/nightly/html/classvtkWedge.html#details) occurs 628 times.
    Cell type [vtkPyramid](https://www.vtk.org/doc/nightly/html/classvtkPyramid.html#details) occurs 4376 times.
The clipped dataset contains a 
  [vtkUnstructuredGrid](https://www.vtk.org/doc/nightly/html/classvtkUnstructuredGrid.html#details) that has 25655 cells
    Cell type [vtkTetra](https://www.vtk.org/doc/nightly/html/classvtkTetra.html#details) occurs 3715 times.
    Cell type [vtkHexahedron](https://www.vtk.org/doc/nightly/html/classvtkHexahedron.html#details) occurs 16984 times.
    Cell type [vtkWedge](https://www.vtk.org/doc/nightly/html/classvtkWedge.html#details) occurs 616 times.
    Cell type [vtkPyramid](https://www.vtk.org/doc/nightly/html/classvtkPyramid.html#details) occurs 4340 times.

Compare these results with ClipUnstructuredGridWithPlane2 (C++), ClipUnstructuredGridWithPlane2 (Python). Also, the resulting vtkUnstructuredGrid's has 1/4 the number of cells.

usage

ClipUnstructuredGridWithPlane treemesh.vtk

thanks

Thanks to Bane Sullivan for sharing the treemesh.vtk unstructured grid dataset.

Other languages

See (Cxx)

Question

If you have a question about this example, please use the VTK Discourse Forum

Code

ClipUnstructuredGridWithPlane.py

#!/usr/bin/env python

import collections

# noinspection PyUnresolvedReferences
import vtkmodules.vtkInteractionStyle
# noinspection PyUnresolvedReferences
import vtkmodules.vtkRenderingOpenGL2
from vtkmodules.vtkCommonColor import vtkNamedColors
from vtkmodules.vtkCommonDataModel import (
    vtkCellTypes,
    vtkPlane
)
from vtkmodules.vtkCommonTransforms import vtkTransform
from vtkmodules.vtkFiltersGeneral import vtkTableBasedClipDataSet
from vtkmodules.vtkIOLegacy import vtkUnstructuredGridReader
from vtkmodules.vtkRenderingCore import (
    vtkActor,
    vtkDataSetMapper,
    vtkRenderWindow,
    vtkRenderWindowInteractor,
    vtkRenderer
)


def get_program_parameters():
    import argparse
    description = 'Use a vtkTableBasedClipDataSet to clip a vtkUnstructuredGrid.'
    epilogue = '''
 Use a vtkTableBasedClipDataSet to clip a vtkUnstructuredGrid.
 The resulting output and clipped output are presented in yellow and red respectively.
 To illustrate the clipped interfaces, the example uses a vtkTransform to rotate each
    output about their centers.
 Note: This clipping filter does retain the original cells if they are not clipped.

   '''
    parser = argparse.ArgumentParser(description=description, epilog=epilogue,
                                     formatter_class=argparse.RawDescriptionHelpFormatter)
    parser.add_argument('filename', help='treemesh.vtk')
    args = parser.parse_args()
    return args.filename


def main():
    filename = get_program_parameters()

    # Create the reader for the data.
    reader = vtkUnstructuredGridReader()
    reader.SetFileName(filename)
    reader.Update()

    bounds = reader.GetOutput().GetBounds()
    center = reader.GetOutput().GetCenter()

    colors = vtkNamedColors()
    renderer = vtkRenderer()
    renderer.SetBackground(colors.GetColor3d('Wheat'))
    renderer.UseHiddenLineRemovalOn()

    renderWindow = vtkRenderWindow()
    renderWindow.AddRenderer(renderer)
    renderWindow.SetSize(640, 480)

    interactor = vtkRenderWindowInteractor()
    interactor.SetRenderWindow(renderWindow)

    xnorm = [-1.0, -1.0, 1.0]

    clipPlane = vtkPlane()
    clipPlane.SetOrigin(reader.GetOutput().GetCenter())
    clipPlane.SetNormal(xnorm)

    clipper = vtkTableBasedClipDataSet()
    clipper.SetClipFunction(clipPlane)
    clipper.SetInputData(reader.GetOutput())
    clipper.SetValue(0.0)
    clipper.GenerateClippedOutputOn()
    clipper.Update()

    insideMapper = vtkDataSetMapper()
    insideMapper.SetInputData(clipper.GetOutput())
    insideMapper.ScalarVisibilityOff()

    insideActor = vtkActor()
    insideActor.SetMapper(insideMapper)
    insideActor.GetProperty().SetDiffuseColor(colors.GetColor3d('banana'))
    insideActor.GetProperty().SetAmbient(.3)
    insideActor.GetProperty().EdgeVisibilityOn()

    clippedMapper = vtkDataSetMapper()
    clippedMapper.SetInputData(clipper.GetClippedOutput())
    clippedMapper.ScalarVisibilityOff()

    clippedActor = vtkActor()
    clippedActor.SetMapper(clippedMapper)
    clippedActor.GetProperty().SetDiffuseColor(colors.GetColor3d('tomato'))
    insideActor.GetProperty().SetAmbient(.3)
    clippedActor.GetProperty().EdgeVisibilityOn()

    # Create transforms to make a better visualization
    insideTransform = vtkTransform()
    insideTransform.Translate(-(bounds[1] - bounds[0]) * 0.75, 0, 0)
    insideTransform.Translate(center[0], center[1], center[2])
    insideTransform.RotateY(-120.0)
    insideTransform.Translate(-center[0], -center[1], -center[2])
    insideActor.SetUserTransform(insideTransform)

    clippedTransform = vtkTransform()
    clippedTransform.Translate((bounds[1] - bounds[0]) * 0.75, 0, 0)
    clippedTransform.Translate(center[0], center[1], center[2])
    clippedTransform.RotateY(60.0)
    clippedTransform.Translate(-center[0], -center[1], -center[2])
    clippedActor.SetUserTransform(clippedTransform)

    renderer.AddViewProp(clippedActor)
    renderer.AddViewProp(insideActor)

    renderer.ResetCamera()
    renderer.GetActiveCamera().Dolly(1.4)
    renderer.ResetCameraClippingRange()
    renderWindow.Render()
    renderWindow.SetWindowName('ClipUnstructuredGridWithPlane')
    renderWindow.Render()

    interactor.Start()

    # Generate a report
    numberOfCells = clipper.GetOutput().GetNumberOfCells()
    print('------------------------')
    print('The inside dataset contains a \n', clipper.GetOutput().GetClassName(), ' that has ', numberOfCells, ' cells')
    cellMap = dict()
    for i in range(0, numberOfCells):
        cellMap.setdefault(clipper.GetOutput().GetCellType(i), 0)
        cellMap[clipper.GetOutput().GetCellType(i)] += 1
    # Sort by key and put into an OrderedDict.
    # An OrderedDict remembers the order in which the keys have been inserted.
    for k, v in collections.OrderedDict(sorted(cellMap.items())).items():
        print('\tCell type ', vtkCellTypes.GetClassNameFromTypeId(k), ' occurs ', v, ' times.')

    numberOfCells = clipper.GetClippedOutput().GetNumberOfCells()
    print('------------------------')
    print('The clipped dataset contains a \n', clipper.GetClippedOutput().GetClassName(), ' that has ', numberOfCells,
          ' cells')
    outsideCellMap = dict()
    for i in range(0, numberOfCells):
        outsideCellMap.setdefault(clipper.GetClippedOutput().GetCellType(i), 0)
        outsideCellMap[clipper.GetClippedOutput().GetCellType(i)] += 1
    for k, v in collections.OrderedDict(sorted(outsideCellMap.items())).items():
        print('\tCell type ', vtkCellTypes.GetClassNameFromTypeId(k), ' occurs ', v, ' times.')


if __name__ == '__main__':
    main()