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
|
#!/usr/bin/env python
## Program: VMTK
## Module: $RCSfile: vmtksurfacecapper.py,v $
## Language: Python
## Date: $Date: 2006/07/17 09:53:14 $
## Version: $Revision: 1.8 $
## 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 sys
import vtk
import vtkvmtk
import pypes
vmtksurfacecapper = 'vmtkSurfaceCapper'
class vmtkSurfaceCapper(pypes.pypeScript):
def __init__(self):
pypes.pypeScript.__init__(self)
self.Surface = None
self.TriangleOutput = 1
self.CellEntityIdsArrayName = 'CellEntityIds'
self.CellEntityIdOffset = 1
self.Interactive = 1
self.Method = 'simple'
self.ConstraintFactor = 1.0
self.NumberOfRings = 8
self.vmtkRenderer = None
self.OwnRenderer = 0
self.SetScriptName('vmtksurfacecapper')
self.SetScriptDoc('add caps to the holes of a surface, assigning an id to each cap for easy specification of boundary conditions ("simple" method only).')
self.SetInputMembers([
['Surface','i','vtkPolyData',1,'','the input surface','vmtksurfacereader'],
['Method','method','str',1,'["simple","centerpoint","smooth","annular"]','capping method'],
['TriangleOutput','triangle','bool',1,'','toggle triangulation of the output'],
['CellEntityIdsArrayName','entityidsarray','str',1,'','name of the array where the id of the caps have to be stored'],
['CellEntityIdOffset','entityidoffset','int',1,'(0,)','offset for entity ids ("simple" method only")'],
['ConstraintFactor','constraint','float',1,'','amount of influence of the shape of the surface near the boundary on the shape of the cap ("smooth" method only)'],
['NumberOfRings','rings','int',1,'(0,)','number of rings composing the cap ("smooth" method only)'],
['Interactive','interactive','bool',1],
['vmtkRenderer','renderer','vmtkRenderer',1,'','external renderer']
])
self.SetOutputMembers([
['Surface','o','vtkPolyData',1,'','the output surface','vmtksurfacewriter'],
['CellEntityIdsArrayName','entityidsarray','str',1,'','name of the array where the id of the caps are stored']
])
def LabelValidator(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 self.Surface == None:
self.PrintError('Error: No input surface.')
# cleaner = vtk.vtkCleanPolyData()
# cleaner.SetInput(self.Surface)
# cleaner.Update()
#
# triangleFilter = vtk.vtkTriangleFilter()
# triangleFilter.SetInput(cleaner.GetOutput())
# triangleFilter.Update()
#
# self.Surface = triangleFilter.GetOutput()
boundaryIds = vtk.vtkIdList()
if self.Interactive:
if not self.vmtkRenderer:
import vmtkrenderer
self.vmtkRenderer = vmtkrenderer.vmtkRenderer()
self.vmtkRenderer.Initialize()
self.OwnRenderer = 1
self.vmtkRenderer.RegisterScript(self)
boundaryExtractor = vtkvmtk.vtkvmtkPolyDataBoundaryExtractor()
boundaryExtractor.SetInput(self.Surface)
boundaryExtractor.Update()
boundaries = boundaryExtractor.GetOutput()
numberOfBoundaries = boundaries.GetNumberOfCells()
seedPoints = vtk.vtkPoints()
for i in range(numberOfBoundaries):
barycenter = [0.0, 0.0, 0.0]
vtkvmtk.vtkvmtkBoundaryReferenceSystems.ComputeBoundaryBarycenter(boundaries.GetCell(i).GetPoints(),barycenter)
seedPoints.InsertNextPoint(barycenter)
seedPolyData = vtk.vtkPolyData()
seedPolyData.SetPoints(seedPoints)
seedPolyData.Update()
labelsMapper = vtk.vtkLabeledDataMapper();
labelsMapper.SetInput(seedPolyData)
labelsMapper.SetLabelModeToLabelIds()
labelsActor = vtk.vtkActor2D()
labelsActor.SetMapper(labelsMapper)
self.vmtkRenderer.Renderer.AddActor(labelsActor)
surfaceMapper = vtk.vtkPolyDataMapper()
surfaceMapper.SetInput(self.Surface)
surfaceMapper.ScalarVisibilityOff()
surfaceActor = vtk.vtkActor()
surfaceActor.SetMapper(surfaceMapper)
surfaceActor.GetProperty().SetOpacity(0.25)
self.vmtkRenderer.Renderer.AddActor(surfaceActor)
#self.vmtkRenderer.Render()
#self.vmtkRenderer.Renderer.RemoveActor(labelsActor)
#self.vmtkRenderer.Renderer.RemoveActor(surfaceActor)
ok = False
while not ok:
labelString = self.InputText("Please input boundary ids: ",self.LabelValidator)
labels = [int(label) for label in labelString.split()]
ok = True
for label in labels:
if label not in range(numberOfBoundaries):
ok = False
for label in labels:
boundaryIds.InsertNextId(label)
if self.Method == 'simple':
capper = vtkvmtk.vtkvmtkSimpleCapPolyData()
capper.SetInput(self.Surface)
if self.Interactive:
capper.SetBoundaryIds(boundaryIds)
capper.SetCellEntityIdsArrayName(self.CellEntityIdsArrayName)
capper.SetCellEntityIdOffset(self.CellEntityIdOffset)
capper.Update()
self.Surface = capper.GetOutput()
elif self.Method == 'centerpoint':
capper = vtkvmtk.vtkvmtkCapPolyData()
capper.SetInput(self.Surface)
if self.Interactive:
capper.SetBoundaryIds(boundaryIds)
capper.SetDisplacement(0.0)
capper.SetInPlaneDisplacement(0.0)
capper.Update()
self.Surface = capper.GetOutput()
elif self.Method == 'smooth':
triangle = vtk.vtkTriangleFilter()
triangle.SetInput(self.Surface)
triangle.PassLinesOff()
triangle.PassVertsOff()
triangle.Update()
capper = vtkvmtk.vtkvmtkSmoothCapPolyData()
capper.SetInput(triangle.GetOutput())
capper.SetConstraintFactor(self.ConstraintFactor)
capper.SetNumberOfRings(self.NumberOfRings)
if self.Interactive:
capper.SetBoundaryIds(boundaryIds)
capper.Update()
self.Surface = capper.GetOutput()
elif self.Method == 'annular':
capper = vtkvmtk.vtkvmtkAnnularCapPolyData()
capper.SetInput(self.Surface)
capper.SetCellEntityIdsArrayName(self.CellEntityIdsArrayName)
capper.SetCellEntityIdOffset(self.CellEntityIdOffset)
capper.Update()
self.Surface = capper.GetOutput()
if self.TriangleOutput == 1:
triangle = vtk.vtkTriangleFilter()
triangle.SetInput(self.Surface)
triangle.PassLinesOff()
triangle.PassVertsOff()
triangle.Update()
self.Surface = triangle.GetOutput()
normals = vtk.vtkPolyDataNormals()
normals.SetInput(self.Surface)
normals.AutoOrientNormalsOn()
normals.SplittingOff()
normals.ConsistencyOn()
normals.Update()
self.Surface = normals.GetOutput()
if self.Surface.GetSource():
self.Surface.GetSource().UnRegisterAllOutputs()
if __name__=='__main__':
main = pypes.pypeMain()
main.Arguments = sys.argv
main.Execute()
|