Skip to content

HTTPS clone URL

Subversion checkout URL

You can clone with HTTPS or Subversion.

Download ZIP
Newer
Older
100644 211 lines (174 sloc) 8.117 kb
665c986 Initial import
lucantiga authored
1 #!/usr/bin/env python
2
3 ## Program: VMTK
4 ## Module: $RCSfile: vmtksurfacecapper.py,v $
5 ## Language: Python
6 ## Date: $Date: 2006/07/17 09:53:14 $
7 ## Version: $Revision: 1.8 $
8
9 ## Copyright (c) Luca Antiga, David Steinman. All rights reserved.
10 ## See LICENCE file for details.
11
12 ## This software is distributed WITHOUT ANY WARRANTY; without even
13 ## the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
14 ## PURPOSE. See the above copyright notices for more information.
15
16
17 import sys
18 import vtk
19
2162ceb Fixed Python module layout and started using CPack
lucantiga authored
20 import vtkvmtk
665c986 Initial import
lucantiga authored
21 import pypes
22
23 vmtksurfacecapper = 'vmtkSurfaceCapper'
24
25 class vmtkSurfaceCapper(pypes.pypeScript):
26
27 def __init__(self):
28
29 pypes.pypeScript.__init__(self)
30
31 self.Surface = None
32 self.TriangleOutput = 1
00517f2 Fixed bugs in surface remeshing that lead to crashes in particular cases
lucantiga authored
33 self.CellEntityIdsArrayName = 'CellEntityIds'
988c9a8 Added vmtkmeshgenerator with related classes. Fixed bug in vmtkimagevois...
lucantiga authored
34 self.CellEntityIdOffset = 1
665c986 Initial import
lucantiga authored
35 self.Interactive = 1
2c08c24 Added smooth capping code
lucantiga authored
36 self.Method = 'simple'
37 self.ConstraintFactor = 1.0
38 self.NumberOfRings = 8
665c986 Initial import
lucantiga authored
39
40 self.vmtkRenderer = None
41 self.OwnRenderer = 0
42
43 self.SetScriptName('vmtksurfacecapper')
2c08c24 Added smooth capping code
lucantiga authored
44 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).')
665c986 Initial import
lucantiga authored
45 self.SetInputMembers([
cbda38c Several improvements to PypeS for better Slicer3 integration and vmtk co...
lucantiga authored
46 ['Surface','i','vtkPolyData',1,'','the input surface','vmtksurfacereader'],
244bca2 Added annular capping. Starting to integrate and generalize Marina's cod...
lucantiga authored
47 ['Method','method','str',1,'["simple","centerpoint","smooth","annular"]','capping method'],
cbda38c Several improvements to PypeS for better Slicer3 integration and vmtk co...
lucantiga authored
48 ['TriangleOutput','triangle','bool',1,'','toggle triangulation of the output'],
00517f2 Fixed bugs in surface remeshing that lead to crashes in particular cases
lucantiga authored
49 ['CellEntityIdsArrayName','entityidsarray','str',1,'','name of the array where the id of the caps have to be stored'],
988c9a8 Added vmtkmeshgenerator with related classes. Fixed bug in vmtkimagevois...
lucantiga authored
50 ['CellEntityIdOffset','entityidoffset','int',1,'(0,)','offset for entity ids ("simple" method only")'],
cbda38c Several improvements to PypeS for better Slicer3 integration and vmtk co...
lucantiga authored
51 ['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)'],
4676628 Added numeric range checking to pypes
lucantiga authored
52 ['NumberOfRings','rings','int',1,'(0,)','number of rings composing the cap ("smooth" method only)'],
cbda38c Several improvements to PypeS for better Slicer3 integration and vmtk co...
lucantiga authored
53 ['Interactive','interactive','bool',1],
54 ['vmtkRenderer','renderer','vmtkRenderer',1,'','external renderer']
665c986 Initial import
lucantiga authored
55 ])
56 self.SetOutputMembers([
9f12639 Added surface remeshing script
lucantiga authored
57 ['Surface','o','vtkPolyData',1,'','the output surface','vmtksurfacewriter'],
00517f2 Fixed bugs in surface remeshing that lead to crashes in particular cases
lucantiga authored
58 ['CellEntityIdsArrayName','entityidsarray','str',1,'','name of the array where the id of the caps are stored']
9f12639 Added surface remeshing script
lucantiga authored
59 ])
665c986 Initial import
lucantiga authored
60
61 def LabelValidator(self,text):
62 import string
63 if not text:
64 return 0
707e538 Fixed flow extensions to respect number of boundary points. Improved val...
lucantiga authored
65 if not text.split():
66 return 0
665c986 Initial import
lucantiga authored
67 for char in text:
68 if char not in string.digits + " ":
69 return 0
70 return 1
71
72 def Execute(self):
73
74 if self.Surface == None:
75 self.PrintError('Error: No input surface.')
76
ba18837 git-svn-id: https://vmtk.svn.sourceforge.net/svnroot/vmtk@81 94be75cf-be...
lucantiga authored
77 # cleaner = vtk.vtkCleanPolyData()
78 # cleaner.SetInput(self.Surface)
79 # cleaner.Update()
80 #
81 # triangleFilter = vtk.vtkTriangleFilter()
82 # triangleFilter.SetInput(cleaner.GetOutput())
83 # triangleFilter.Update()
84 #
85 # self.Surface = triangleFilter.GetOutput()
86
665c986 Initial import
lucantiga authored
87 boundaryIds = vtk.vtkIdList()
88
89 if self.Interactive:
90 if not self.vmtkRenderer:
91 import vmtkrenderer
92 self.vmtkRenderer = vmtkrenderer.vmtkRenderer()
93 self.vmtkRenderer.Initialize()
94 self.OwnRenderer = 1
95
82a46ff @SaraZanchi Modified Scripts
SaraZanchi authored
96 self.vmtkRenderer.RegisterScript(self)
97
665c986 Initial import
lucantiga authored
98 boundaryExtractor = vtkvmtk.vtkvmtkPolyDataBoundaryExtractor()
99 boundaryExtractor.SetInput(self.Surface)
100 boundaryExtractor.Update()
1d8aaa2 Added flow extensions based on normal to boundary and specified radius, ...
lucantiga authored
101
665c986 Initial import
lucantiga authored
102 boundaries = boundaryExtractor.GetOutput()
103 numberOfBoundaries = boundaries.GetNumberOfCells()
104 seedPoints = vtk.vtkPoints()
105 for i in range(numberOfBoundaries):
106 barycenter = [0.0, 0.0, 0.0]
107 vtkvmtk.vtkvmtkBoundaryReferenceSystems.ComputeBoundaryBarycenter(boundaries.GetCell(i).GetPoints(),barycenter)
108 seedPoints.InsertNextPoint(barycenter)
109 seedPolyData = vtk.vtkPolyData()
110 seedPolyData.SetPoints(seedPoints)
111 seedPolyData.Update()
112 labelsMapper = vtk.vtkLabeledDataMapper();
113 labelsMapper.SetInput(seedPolyData)
114 labelsMapper.SetLabelModeToLabelIds()
115 labelsActor = vtk.vtkActor2D()
116 labelsActor.SetMapper(labelsMapper)
117
118 self.vmtkRenderer.Renderer.AddActor(labelsActor)
119
120 surfaceMapper = vtk.vtkPolyDataMapper()
121 surfaceMapper.SetInput(self.Surface)
122 surfaceMapper.ScalarVisibilityOff()
123 surfaceActor = vtk.vtkActor()
124 surfaceActor.SetMapper(surfaceMapper)
125 surfaceActor.GetProperty().SetOpacity(0.25)
126
127 self.vmtkRenderer.Renderer.AddActor(surfaceActor)
128
82a46ff @SaraZanchi Modified Scripts
SaraZanchi authored
129 #self.vmtkRenderer.Render()
130
131 #self.vmtkRenderer.Renderer.RemoveActor(labelsActor)
132 #self.vmtkRenderer.Renderer.RemoveActor(surfaceActor)
133
665c986 Initial import
lucantiga authored
134 ok = False
135 while not ok:
136 labelString = self.InputText("Please input boundary ids: ",self.LabelValidator)
137 labels = [int(label) for label in labelString.split()]
138 ok = True
139 for label in labels:
140 if label not in range(numberOfBoundaries):
141 ok = False
142
143 for label in labels:
144 boundaryIds.InsertNextId(label)
145
2c08c24 Added smooth capping code
lucantiga authored
146 if self.Method == 'simple':
147 capper = vtkvmtk.vtkvmtkSimpleCapPolyData()
148 capper.SetInput(self.Surface)
149 if self.Interactive:
150 capper.SetBoundaryIds(boundaryIds)
9f12639 Added surface remeshing script
lucantiga authored
151 capper.SetCellEntityIdsArrayName(self.CellEntityIdsArrayName)
988c9a8 Added vmtkmeshgenerator with related classes. Fixed bug in vmtkimagevois...
lucantiga authored
152 capper.SetCellEntityIdOffset(self.CellEntityIdOffset)
2c08c24 Added smooth capping code
lucantiga authored
153 capper.Update()
154 self.Surface = capper.GetOutput()
155 elif self.Method == 'centerpoint':
156 capper = vtkvmtk.vtkvmtkCapPolyData()
157 capper.SetInput(self.Surface)
158 if self.Interactive:
159 capper.SetBoundaryIds(boundaryIds)
160 capper.SetDisplacement(0.0)
161 capper.SetInPlaneDisplacement(0.0)
162 capper.Update()
163 self.Surface = capper.GetOutput()
164 elif self.Method == 'smooth':
165 triangle = vtk.vtkTriangleFilter()
166 triangle.SetInput(self.Surface)
167 triangle.PassLinesOff()
168 triangle.PassVertsOff()
169 triangle.Update()
170 capper = vtkvmtk.vtkvmtkSmoothCapPolyData()
171 capper.SetInput(triangle.GetOutput())
172 capper.SetConstraintFactor(self.ConstraintFactor)
173 capper.SetNumberOfRings(self.NumberOfRings)
174 if self.Interactive:
175 capper.SetBoundaryIds(boundaryIds)
176 capper.Update()
177 self.Surface = capper.GetOutput()
244bca2 Added annular capping. Starting to integrate and generalize Marina's cod...
lucantiga authored
178 elif self.Method == 'annular':
179 capper = vtkvmtk.vtkvmtkAnnularCapPolyData()
180 capper.SetInput(self.Surface)
181 capper.SetCellEntityIdsArrayName(self.CellEntityIdsArrayName)
182 capper.SetCellEntityIdOffset(self.CellEntityIdOffset)
183 capper.Update()
184 self.Surface = capper.GetOutput()
665c986 Initial import
lucantiga authored
185
186 if self.TriangleOutput == 1:
187 triangle = vtk.vtkTriangleFilter()
188 triangle.SetInput(self.Surface)
189 triangle.PassLinesOff()
190 triangle.PassVertsOff()
191 triangle.Update()
192 self.Surface = triangle.GetOutput()
193
194 normals = vtk.vtkPolyDataNormals()
195 normals.SetInput(self.Surface)
196 normals.AutoOrientNormalsOn()
197 normals.SplittingOff()
198 normals.ConsistencyOn()
199 normals.Update()
200 self.Surface = normals.GetOutput()
201
202 if self.Surface.GetSource():
203 self.Surface.GetSource().UnRegisterAllOutputs()
204
205
206 if __name__=='__main__':
207
208 main = pypes.pypeMain()
209 main.Arguments = sys.argv
210 main.Execute()
Something went wrong with that request. Please try again.