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
|
#!/usr/bin/env python
## Program: VMTK
## Module: $RCSfile: vmtkmeshbranchclipper.py,v $
## Language: Python
## Date: $Date: 2006/02/23 09:31:39 $
## Version: $Revision: 1.1 $
## 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 vtkvmtk
import sys
import pypes
import vmtkrenderer
import vmtkcenterlineviewer
vmtkmeshbranchclipper = 'vmtkMeshBranchClipper'
class vmtkMeshBranchClipper(pypes.pypeScript):
def __init__(self):
pypes.pypeScript.__init__(self)
self.Mesh = None
self.Centerlines = None
self.RadiusArrayName = ''
self.CutoffRadiusFactor = 1E16
self.ClipValue = 0.0
self.BlankingArrayName = ''
self.GroupIdsArrayName = ''
self.GroupIds = []
self.InsideOut = 0
self.UseRadiusInformation = 1
self.vmtkRenderer = None
self.OwnRenderer = 0
self.Interactive = 0
self.SetScriptName('vmtkmeshbranchclipper')
self.SetInputMembers([
['Mesh','i','vtkUnstructuredGrid',1,'','','vmtkmeshreader'],
['Centerlines','centerlines','vtkPolyData',1,'','','vmtksurfacereader'],
['GroupIdsArrayName','groupidsarray','str',1],
['GroupIds','groupids','int',-1],
['InsideOut','insideout','bool',1],
['UseRadiusInformation','useradius','bool',1],
['RadiusArrayName','radiusarray','str',1],
['BlankingArrayName','blankingarray','str',1],
['CutoffRadiusFactor','cutoffradiusfactor','float',1,'(0.0,)'],
['ClipValue','clipvalue','float',1],
['Interactive','interactive','bool',1],
['vmtkRenderer','renderer','vmtkRenderer',1,'','external renderer']
])
self.SetOutputMembers([
['Mesh','o','vtkUnstructuredGrid',1,'','','vmtkmeshwriter'],
['Centerlines','ocenterlines','vtkPolyData',1,'','','vmtksurfacewriter']
])
def GroupIdsValidator(self,text):
import string
if not text:
return 0
if not text.split():
return 0
for char in text:
if char not in string.digits + " ":
return 0
return 1
def Execute(self):
if not self.Mesh:
self.PrintError('Error: No input mesh.')
if not self.Centerlines:
self.PrintError('Error: No input centerlines.')
if self.Interactive and not self.vmtkRenderer:
self.vmtkRenderer = vmtkrenderer.vmtkRenderer()
self.vmtkRenderer.Initialize()
self.OwnRenderer = 1
if self.Interactive:
self.vmtkRenderer.RegisterScript(self)
viewer = vmtkcenterlineviewer.vmtkCenterlineViewer()
viewer.Centerlines = self.Centerlines
viewer.CellDataArrayName = self.GroupIdsArrayName
viewer.vmtkRenderer = self.vmtkRenderer
viewer.InputText = self.InputText
viewer.OutputText = self.OutputText
viewer.PrintError = self.PrintError
viewer.PringLog = self.PrintLog
viewer.Display = 0
viewer.Execute()
groupIdsString = self.InputText("Please input groupIds to clip:\n",self.GroupIdsValidator)
self.GroupIds = [int(groupId) for groupId in groupIdsString.split()]
clipper = vtkvmtk.vtkvmtkUnstructuredGridCenterlineGroupsClipper()
clipper.SetInput(self.Mesh)
clipper.SetCenterlines(self.Centerlines)
clipper.SetCenterlineGroupIdsArrayName(self.GroupIdsArrayName)
clipper.SetGroupIdsArrayName(self.GroupIdsArrayName)
clipper.SetCenterlineRadiusArrayName(self.RadiusArrayName)
clipper.SetBlankingArrayName(self.BlankingArrayName)
clipper.SetCutoffRadiusFactor(self.CutoffRadiusFactor)
clipper.SetClipValue(self.ClipValue)
clipper.SetUseRadiusInformation(self.UseRadiusInformation)
if self.GroupIds:
groupIds = vtk.vtkIdList()
for groupId in self.GroupIds:
groupIds.InsertNextId(groupId)
clipper.SetCenterlineGroupIds(groupIds)
clipper.ClipAllCenterlineGroupIdsOff()
else:
clipper.ClipAllCenterlineGroupIdsOn()
if not self.InsideOut:
clipper.GenerateClippedOutputOff()
else:
clipper.GenerateClippedOutputOn()
clipper.Update()
if not self.InsideOut:
self.Mesh = clipper.GetOutput()
else:
self.Mesh = clipper.GetClippedOutput()
if self.Mesh:
if self.Mesh.GetSource():
self.Mesh.GetSource().UnRegisterAllOutputs()
if self.Centerlines.GetSource():
self.Centerlines.GetSource().UnRegisterAllOutputs()
if self.OwnRenderer:
self.vmtkRenderer.Deallocate()
if __name__=='__main__':
main = pypes.pypeMain()
main.Arguments = sys.argv
main.Execute()
|