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 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291
|
#!/usr/bin/env python
## Program: VMTK
## Module: $RCSfile: vmtksurfaceendclipper.py,v $
## Language: Python
## Date: $Date: 2015/12/01 14:45:13 $
## Version: $Revision: 1.0 $
## 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.
## Note: this class was contributed by
## Simone Manini
## Orobix Srl
import vtk
import sys
import vtkvmtk
import vmtkrenderer
import pypes
vmtksurfaceendclipper = 'vmtkSurfaceEndClipper'
## TODO: make SeedSelector a separate pype script to be used in other contexts
class vmtkSeedSelector(object):
def __init__(self):
self._Surface = None
self._SeedIds = None
self._SourceSeedIds = vtk.vtkIdList()
self._TargetSeedIds = vtk.vtkIdList()
self.PrintError = None
self.PrintLog = None
self.InputText = None
self.OutputText = None
self.InputInfo = None
def SetSurface(self,surface):
self._Surface = surface
def GetSurface(self):
return self._Surface
def GetSourceSeedIds(self):
return self._SourceSeedIds
def GetTargetSeedIds(self):
return self._TargetSeedIds
def Execute(self):
pass
class vmtkPickPointSeedSelector(vmtkSeedSelector):
def __init__(self):
vmtkSeedSelector.__init__(self)
self.PickedSeedIds = vtk.vtkIdList()
self.PickedSeeds = vtk.vtkPolyData()
self.vmtkRenderer = None
self.OwnRenderer = 0
self.Script = None
def UndoCallback(self, obj):
self.InitializeSeeds()
self.PickedSeeds.Modified()
self.vmtkRenderer.RenderWindow.Render()
def PickCallback(self, obj):
picker = vtk.vtkCellPicker()
picker.SetTolerance(1E-4 * self._Surface.GetLength())
eventPosition = self.vmtkRenderer.RenderWindowInteractor.GetEventPosition()
result = picker.Pick(float(eventPosition[0]),float(eventPosition[1]),0.0,self.vmtkRenderer.Renderer)
if result == 0:
return
pickPosition = picker.GetPickPosition()
pickedCellPointIds = self._Surface.GetCell(picker.GetCellId()).GetPointIds()
minDistance = 1E10
pickedSeedId = -1
for i in range(pickedCellPointIds.GetNumberOfIds()):
distance = vtk.vtkMath.Distance2BetweenPoints(pickPosition,self._Surface.GetPoint(pickedCellPointIds.GetId(i)))
if distance < minDistance:
minDistance = distance
pickedSeedId = pickedCellPointIds.GetId(i)
if pickedSeedId == -1:
pickedSeedId = pickedCellPointIds.GetId(0)
self.PickedSeedIds.InsertNextId(pickedSeedId)
point = self._Surface.GetPoint(pickedSeedId)
self.PickedSeeds.GetPoints().InsertNextPoint(point)
self.PickedSeeds.Modified()
self.vmtkRenderer.RenderWindow.Render()
def InitializeSeeds(self):
self.PickedSeedIds.Initialize()
self.PickedSeeds.Initialize()
seedPoints = vtk.vtkPoints()
self.PickedSeeds.SetPoints(seedPoints)
def Execute(self):
if (self._Surface == None):
self.PrintError('vmtkPickPointSeedSelector Error: Surface not set.')
return
self._SourceSeedIds.Initialize()
self._TargetSeedIds.Initialize()
if not self.vmtkRenderer:
self.vmtkRenderer = vmtkrenderer.vmtkRenderer()
self.vmtkRenderer.Initialize()
self.OwnRenderer = 1
self.vmtkRenderer.RegisterScript(self.Script)
glyphs = vtk.vtkGlyph3D()
glyphSource = vtk.vtkSphereSource()
glyphs.SetInputData(self.PickedSeeds)
glyphs.SetSourceConnection(glyphSource.GetOutputPort())
glyphs.SetScaleModeToDataScalingOff()
glyphs.SetScaleFactor(self._Surface.GetLength()*0.01)
glyphMapper = vtk.vtkPolyDataMapper()
glyphMapper.SetInputConnection(glyphs.GetOutputPort())
self.SeedActor = vtk.vtkActor()
self.SeedActor.SetMapper(glyphMapper)
self.SeedActor.GetProperty().SetColor(1.0,0.0,0.0)
self.SeedActor.PickableOff()
self.vmtkRenderer.Renderer.AddActor(self.SeedActor)
##self.vmtkRenderer.RenderWindowInteractor.AddObserver("KeyPressEvent", self.KeyPressed)
self.vmtkRenderer.AddKeyBinding('u','Undo.',self.UndoCallback)
self.vmtkRenderer.AddKeyBinding('space','Add points.',self.PickCallback)
surfaceMapper = vtk.vtkPolyDataMapper()
surfaceMapper.SetInputData(self._Surface)
surfaceMapper.ScalarVisibilityOff()
surfaceActor = vtk.vtkActor()
surfaceActor.SetMapper(surfaceMapper)
surfaceActor.GetProperty().SetOpacity(1.0)
self.vmtkRenderer.Renderer.AddActor(surfaceActor)
self.InputInfo('Please position the mouse and press space to add a point on the main body of the model, \'u\' to undo\n')
any = 0
while any == 0:
self.InitializeSeeds()
self.vmtkRenderer.Render()
any = self.PickedSeedIds.GetNumberOfIds()
self._SourceSeedIds.DeepCopy(self.PickedSeedIds)
self.InputInfo('Please position the mouse and press space to add points at clip locations, \'u\' to undo\n')
any = 0
while any == 0:
self.InitializeSeeds()
self.vmtkRenderer.Render()
any = self.PickedSeedIds.GetNumberOfIds()
self._TargetSeedIds.DeepCopy(self.PickedSeedIds)
if self.OwnRenderer:
self.vmtkRenderer.Deallocate()
class vmtkSurfaceEndClipper(pypes.pypeScript):
def __init__(self):
pypes.pypeScript.__init__(self)
self.Surface = None
self.vmtkRenderer = None
self.OwnRenderer = 0
self.SeedSelector = None
self.SeedSelectorName = 'pointlist'
self.SourcePoints = []
self.TargetPoints = []
self.Actor = None
self.CleanOutput = 1
self.SetScriptName('vmtksurfaceendclipper')
self.SetScriptDoc('interactively clip a tubular surface with normals estimated at seed locations')
self.SetInputMembers([
['Surface','i','vtkPolyData',1,'','the input surface','vmtksurfacereader'],
['vmtkRenderer','renderer','vmtkRenderer',1,'','external renderer']
])
self.SetOutputMembers([
['Surface','o','vtkPolyData',1,'','the output surface','vmtksurfacewriter'],
])
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.SeedSelector = vmtkPickPointSeedSelector()
self.SeedSelector.vmtkRenderer = self.vmtkRenderer
self.SeedSelector.Script = self
self.SeedSelector.SetSurface(self.Surface)
self.SeedSelector.InputInfo = self.InputInfo
self.SeedSelector.InputText = self.InputText
self.SeedSelector.OutputText = self.OutputText
self.SeedSelector.PrintError = self.PrintError
self.SeedSelector.PrintLog = self.PrintLog
self.SeedSelector.Execute()
mainBodySeedIds = self.SeedSelector.GetSourceSeedIds()
mainBodySeedId = mainBodySeedIds.GetId(0)
mainBodyPoint = self.Surface.GetPoint(mainBodySeedId)
clipSeedIds = self.SeedSelector.GetTargetSeedIds()
surfaceCleaner = vtk.vtkCleanPolyData()
surfaceCleaner.SetInputData(self.Surface)
surfaceCleaner.Update()
surfaceTriangulator = vtk.vtkTriangleFilter()
surfaceTriangulator.SetInputConnection(surfaceCleaner.GetOutputPort())
surfaceTriangulator.PassLinesOff()
surfaceTriangulator.PassVertsOff()
surfaceTriangulator.Update()
clippedSurface = surfaceTriangulator.GetOutput()
for i in range(clipSeedIds.GetNumberOfIds()):
seedId = clipSeedIds.GetId(i)
locator = vtk.vtkPointLocator()
locator.SetDataSet(clippedSurface)
locator.BuildLocator()
seedPoint = self.Surface.GetPoint(seedId)
seedPointId = locator.FindClosestPoint(seedPoint)
planeEstimator = vtkvmtk.vtkvmtkPolyDataNormalPlaneEstimator()
planeEstimator.SetInputData(clippedSurface)
planeEstimator.SetOriginPointId(seedPointId)
planeEstimator.Update()
plane = vtk.vtkPlane()
plane.SetOrigin(planeEstimator.GetOrigin())
plane.SetNormal(planeEstimator.GetNormal())
seamFilter = vtkvmtk.vtkvmtkTopologicalSeamFilter()
seamFilter.SetInputData(clippedSurface)
seamFilter.SetClosestPoint(seedPoint)
seamFilter.SetSeamScalarsArrayName("SeamScalars")
seamFilter.SetSeamFunction(plane)
clipper = vtk.vtkClipPolyData()
clipper.SetInputConnection(seamFilter.GetOutputPort())
clipper.GenerateClipScalarsOff()
clipper.GenerateClippedOutputOn()
connectivity = vtk.vtkPolyDataConnectivityFilter()
connectivity.SetInputConnection(clipper.GetOutputPort())
connectivity.SetExtractionModeToClosestPointRegion()
connectivity.SetClosestPoint(mainBodyPoint)
surfaceCleaner = vtk.vtkCleanPolyData()
surfaceCleaner.SetInputConnection(connectivity.GetOutputPort())
surfaceCleaner.Update()
surfaceTriangulator = vtk.vtkTriangleFilter()
surfaceTriangulator.SetInputConnection(surfaceCleaner.GetOutputPort())
surfaceTriangulator.PassLinesOff()
surfaceTriangulator.PassVertsOff()
surfaceTriangulator.Update()
#clippedSurface = vtk.vtkPolyData()
#clippedSurface.DeepCopy(surfaceTriangulator.GetOutput())
clippedSurface = surfaceTriangulator.GetOutput()
self.Surface = clippedSurface
if __name__=='__main__':
main = pypes.pypeMain()
main.Arguments = sys.argv
main.Execute()
|