1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176
|
#!/usr/bin/env python
## Program: VMTK
## Module: $RCSfile: vmtksurfaceregiondrawing.py,v $
## Language: Python
## Date: $Date: 2006/05/26 12:35:13 $
## Version: $Revision: 1.9 $
## Copyright (c) Luca Antiga, David Steinman. All rights reserved.
## See LICENCE file for details.
## This software is distributed WITHOUT ANY WARRANTY; without even
## the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
## PURPOSE. See the above copyright notices for more information.
import vtk
import sys
import vmtkrenderer
import pypes
vmtksurfaceregiondrawing = 'vmtkSurfaceRegionDrawing'
class vmtkSurfaceRegionDrawing(pypes.pypeScript):
def __init__(self):
pypes.pypeScript.__init__(self)
self.Surface = None
self.vmtkRenderer = None
self.OwnRenderer = 0
self.Actor = None
self.ContourWidget = None
self.Interpolator = None
self.Binary = 1
self.OutsideValue = 1.0
self.InsideValue = 0.0
self.ContourScalarsArrayName = 'ContourScalars'
self.SetScriptName('vmtksurfaceregiondrawing')
self.SetScriptDoc('draw a closed contour on a surface and generate a distance field on the surface')
self.SetInputMembers([
['Surface','i','vtkPolyData',1,'','the input surface','vmtksurfacereader'],
['Binary','binary','bool',1,'','fill contour with inside value instead of distance to contour'],
['InsideValue','inside','float',1,'','value with which the surface within the contour is filled'],
['OutsideValue','outside','float',1,'','value with which the surface outside the contour is filled'],
['ContourScalarsArrayName','array','str',1,'','the name of the array where the generated scalars are stored'],
['vmtkRenderer','renderer','vmtkRenderer',1,'','external renderer']
])
self.SetOutputMembers([
['Surface','o','vtkPolyData',1,'','the output surface','vmtksurfacewriter']
])
def ScalarsCallback(self, obj):
rep = vtk.vtkOrientedGlyphContourRepresentation.SafeDownCast(self.ContourWidget.GetRepresentation())
pointIds = vtk.vtkIdList()
self.Interpolator.GetContourPointIds(rep,pointIds)
points = vtk.vtkPoints()
points.SetNumberOfPoints(pointIds.GetNumberOfIds())
for i in range(pointIds.GetNumberOfIds()):
pointId = pointIds.GetId(i)
point = self.Surface.GetPoint(pointId)
points.SetPoint(i,point)
selectionFilter = vtk.vtkSelectPolyData()
selectionFilter.SetInputData(self.Surface)
selectionFilter.SetLoop(points)
selectionFilter.GenerateSelectionScalarsOn()
selectionFilter.SetSelectionModeToSmallestRegion()
selectionFilter.Update()
selectionScalars = selectionFilter.GetOutput().GetPointData().GetScalars()
contourScalars = self.Surface.GetPointData().GetArray(self.ContourScalarsArrayName)
for i in range(contourScalars.GetNumberOfTuples()):
selectionValue = selectionScalars.GetTuple1(i)
if self.Binary:
if selectionValue < 0.0:
contourScalars.SetTuple1(i,self.InsideValue)
else:
contourValue = contourScalars.GetTuple1(i)
if (not contourValue < 0.0 and selectionValue < 0.0) or (contourValue < 0.0 and selectionValue < contourValue):
contourScalars.SetTuple1(i,selectionValue)
self.Actor.GetMapper().SetScalarRange(contourScalars.GetRange(0))
self.Surface.Modified()
self.ContourWidget.Initialize()
def DeleteContourCallback(self, obj):
self.ContourWidget.Initialize()
def InteractCallback(self, obj):
if self.ContourWidget.GetEnabled() == 1:
self.ContourWidget.SetEnabled(0)
else:
self.ContourWidget.SetEnabled(1)
def Display(self):
self.vmtkRenderer.Render()
def Execute(self):
if self.Surface == None:
self.PrintError('Error: no Surface.')
if not self.vmtkRenderer:
self.vmtkRenderer = vmtkrenderer.vmtkRenderer()
self.vmtkRenderer.Initialize()
self.OwnRenderer = 1
self.vmtkRenderer.RegisterScript(self)
triangleFilter = vtk.vtkTriangleFilter()
triangleFilter.SetInputData(self.Surface)
triangleFilter.Update()
self.Surface = triangleFilter.GetOutput()
contourScalars = vtk.vtkDoubleArray()
contourScalars.SetNumberOfComponents(1)
contourScalars.SetNumberOfTuples(self.Surface.GetNumberOfPoints())
contourScalars.SetName(self.ContourScalarsArrayName)
contourScalars.FillComponent(0,self.OutsideValue)
self.Surface.GetPointData().AddArray(contourScalars)
self.Surface.GetPointData().SetActiveScalars(self.ContourScalarsArrayName)
mapper = vtk.vtkPolyDataMapper()
mapper.SetInputData(self.Surface)
mapper.ScalarVisibilityOn()
self.Actor = vtk.vtkActor()
self.Actor.SetMapper(mapper)
self.Actor.GetMapper().SetScalarRange(-1.0,0.0)
self.vmtkRenderer.Renderer.AddActor(self.Actor)
self.ContourWidget = vtk.vtkContourWidget()
self.ContourWidget.SetInteractor(self.vmtkRenderer.RenderWindowInteractor)
rep = vtk.vtkOrientedGlyphContourRepresentation.SafeDownCast(self.ContourWidget.GetRepresentation())
rep.GetLinesProperty().SetColor(1, 0.2, 0)
rep.GetLinesProperty().SetLineWidth(3.0)
pointPlacer = vtk.vtkPolygonalSurfacePointPlacer()
pointPlacer.AddProp(self.Actor)
pointPlacer.GetPolys().AddItem(self.Surface)
rep.SetPointPlacer(pointPlacer)
self.Interpolator = vtk.vtkPolygonalSurfaceContourLineInterpolator()
self.Interpolator.GetPolys().AddItem(self.Surface)
rep.SetLineInterpolator(self.Interpolator)
self.vmtkRenderer.AddKeyBinding('space','Generate scalars',self.ScalarsCallback)
self.vmtkRenderer.AddKeyBinding('d','Delete contour',self.DeleteContourCallback)
self.vmtkRenderer.AddKeyBinding('i','Start interaction',self.InteractCallback)
self.Display()
if self.OwnRenderer:
self.vmtkRenderer.Deallocate()
if __name__=='__main__':
main = pypes.pypeMain()
main.Arguments = sys.argv
main.Execute()
|