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 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356
|
#!/usr/bin/env python
## Program: VMTK
## Module: $RCSfile: vmtkgeodesicsurfaceresolution.py,v $
## Language: Python
## Date: $$
## Version: $$
## 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
## Tangui Morvan
## Kalkulo AS
## Simula Research Laboratory
## This allows the user to specify a edge-length array to be used to specify resolution for surface remeshing
## The array is produced by RBF interpolation of values specified by the user by positioning spheres
## This version use the geodesic distance along the surface instead of hte 3D euclidean distance for hte RBF
import vtk
import sys
import vtkvmtk
import vtkvmtkcontrib
import vmtkrenderer
import pypes
vmtkgeodesicsurfaceresolution = 'vmtkGeodesicSurfaceResolution'
class vmtkGeodesicSurfaceResolution(pypes.pypeScript):
def __init__(self):
pypes.pypeScript.__init__(self)
self.Surface = None
self.ResolutionArrayName = 'ResolutionArray'
self.RBFType = 'biharmonic'
self.Spheres = vtk.vtkPolyData()
self.SphereIds = vtk.vtkIdList()
self.vmtkRenderer = None
self.OwnRenderer = 0
self.DisplayArray = False
self.SurfaceMapper = None
self.CurrentSphereId = -1
self.SphereWidget = None
self.Opacity = 1.
self.SpheresActor = None
self.ScalarBarActor = None
self.InteractionMode = 0
self.ExamineSurface = None
self.ExamineSpheres = vtk.vtkPolyData()
self.ExamineSpheresActor = None
self.ExamineText = None
self.SetScriptName('vtksurfaceresolution')
self.SetInputMembers([
['Surface','i','vtkPolyData',1,'','the input surface','vmtksurfacereader'],
['ResolutionArrayName','resolutionarray','str',1,'','array storing the desired edge length'],
['RBFType','rbftype','str',1,'["thinplatespline","biharmonic","triharmonic"]','the type of RBF interpolation'],
['Opacity','opacity','float',1,'(0.0,1.0)','object opacities in the scene'],
['vmtkRenderer','renderer','vmtkRenderer',1,'','external renderer']
])
self.SetOutputMembers([
['Surface','o','vtkPolyData',1,'','','vmtksurfacewriter']
])
def ComputeArray(self):
rbf = vtkvmtkcontrib.vtkvmtkPolyDataGeodesicRBFInterpolation()
rbf.SetSeedIds(self.SphereIds)
seedValues = self.Spheres.GetPointData().GetScalars()
rbf.SetSeedValues(seedValues)
if self.RBFType == "thinplatespline":
rbf.SetRBFTypeToThinPlateSpline()
elif self.RBFType == "biharmonic":
rbf.SetRBFTypeToBiharmonic()
elif self.RBFType == "triharmonic":
rbf.SetRBFTypeToTriharmonic()
rbf.SetInput(self.Surface)
rbf.SetInterpolatedArrayName(self.ResolutionArrayName)
rbf.Update()
return rbf.GetOutput()
def InitializeSpheres(self):
if (self.InteractionMode==0):
self.SphereIds.Initialize()
self.Spheres.Initialize()
seedPoints = vtk.vtkPoints()
self.Spheres.SetPoints(seedPoints)
self.Spheres.GetPointData().Initialize()
seedRadii = vtk.vtkDoubleArray()
self.Spheres.GetPointData().SetScalars(seedRadii)
self.CurrentSphereId = -1
self.SphereWidget.Off()
else:
self.ExamineSpheres.Initialize()
spherePoints = vtk.vtkPoints()
self.ExamineSpheres.SetPoints(spherePoints)
self.ExamineSpheres.GetPointData().Initialize()
sphereRadii = vtk.vtkDoubleArray()
self.ExamineSpheres.GetPointData().SetScalars(sphereRadii)
def PlaceSphere(self):
if self.CurrentSphereId == -1:
return
self.SphereWidget.SetCenter(self.Spheres.GetPoint(self.CurrentSphereId))
self.SphereWidget.SetRadius(self.Spheres.GetPointData().GetScalars().GetValue(self.CurrentSphereId))
def SphereCallback(self,widget,event_string):
if self.CurrentSphereId == -1:
return
minRadius = self.Surface.GetLength()*0.001
if self.SphereWidget.GetRadius() < minRadius:
self.SphereWidget.SetRadius(minRadius)
self.Spheres.GetPointData().GetScalars().SetValue(self.CurrentSphereId,self.SphereWidget.GetRadius())
self.Spheres.Modified()
def UndoCallback(self, obj):
self.InitializeSpheres()
self.Spheres.Modified()
self.vmtkRenderer.RenderWindow.Render()
def PlusCallback(self, obj):
if self.CurrentSphereId != -1:
#increase sphere radius
newval = self.Spheres.GetPointData().GetScalars().GetValue(self.CurrentSphereId) + self.Surface.GetLength()*0.01
self.Spheres.GetPointData().GetScalars().SetValue(self.CurrentSphereId,newval)
self.Spheres.Modified()
self.PlaceSphere()
self.vmtkRenderer.RenderWindow.Render()
def MinusCallback(self, obj):
if self.CurrentSphereId != -1:
#decrease sphere radius
newval = self.Spheres.GetPointData().GetScalars().GetValue(self.CurrentSphereId) - self.Surface.GetLength()*0.01
if newval> 0:
self.Spheres.GetPointData().GetScalars().SetValue(self.CurrentSphereId,newval)
self.Spheres.Modified()
self.PlaceSphere()
self.vmtkRenderer.RenderWindow.Render()
def NextSphereCallback(self, obj):
if self.CurrentSphereId != -1:
self.CurrentSphereId = (self.CurrentSphereId + 1) % self.Spheres.GetNumberOfPoints();
self.PlaceSphere()
self.vmtkRenderer.RenderWindow.Render()
def PreviousSphereCallback(self, obj):
if self.CurrentSphereId != -1:
self.CurrentSphereId = (self.CurrentSphereId - 1) % self.Spheres.GetNumberOfPoints();
self.PlaceSphere()
self.vmtkRenderer.RenderWindow.Render()
def DisplayCallback(self, obj):
self.DisplayArray = not self.DisplayArray
if self.DisplayArray:
self.ExamineSurface = self.ComputeArray()
self.SurfaceMapper.SetInput(self.ExamineSurface)
self.ExamineSurface.GetPointData().SetActiveScalars(self.ResolutionArrayName)
array = self.ExamineSurface.GetPointData().GetScalars()
if (array):
array.Modified()
self.SurfaceMapper.SetScalarRange(array.GetRange(0))
self.ScalarBarActor.VisibilityOn()
else:
self.SurfaceMapper.SetInput(self.Surface)
self.ScalarBarActor.VisibilityOff()
self.SurfaceMapper.SetScalarVisibility(self.DisplayArray)
self.vmtkRenderer.RenderWindow.Render()
def ExmienModeCallback(self, obj):
#Switch beetween examien and interact mode
if self.InteractionMode == 0:
self.InteractionMode = 1
self.ExamineSurface = self.ComputeArray()
#self.SpheresActor.VisibilityOff()
self.SphereWidget.Off()
self.ExamineSpheresActor.VisibilityOn()
self.ExamineText.VisibilityOn()
self.InitializeSpheres()
else:
self.InteractionMode = 0
#Compute the distances
self.SpheresActor.VisibilityOn()
self.ExamineSpheresActor.VisibilityOff()
self.ExamineText.VisibilityOff()
if (self.CurrentSphereId!=-1):
self.SphereWidget.On()
self.vmtkRenderer.RenderWindow.Render()
def PickCallback(self, obj):
picker = vtk.vtkCellPicker()
picker.SetTolerance(1E-4 * self.Surface.GetLength())
eventPosition = self.vmtkRenderer.RenderWindowInteractor.GetEventPosition()
#eventPosition = obj.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
pickedPointId = -1
for i in range(pickedCellPointIds.GetNumberOfIds()):
distance = vtk.vtkMath.Distance2BetweenPoints(pickPosition,self.Surface.GetPoint(pickedCellPointIds.GetId(i)))
if distance < minDistance:
minDistance = distance
pickedPointId = pickedCellPointIds.GetId(i)
if pickedPointId == -1:
pickedPointId = pickedCellPointIds.GetId(0)
pickedPoint = self.Surface.GetPoint(pickedPointId)
if (self.InteractionMode==0):
self.SphereIds.InsertNextId(pickedPointId)
self.CurrentSphereId = self.Spheres.GetPoints().InsertNextPoint(pickedPoint)
self.Spheres.GetPointData().GetScalars().InsertNextValue(self.Surface.GetLength()*0.01)
self.Spheres.Modified()
self.PlaceSphere()
self.SphereWidget.On()
else:
self.ExamineSpheres.GetPoints().InsertNextPoint(pickedPoint)
length = 0.
array = self.ExamineSurface.GetPointData().GetArray(self.ResolutionArrayName)
if (array):
length = array.GetComponent(pickedPointId,0)
self.ExamineSpheres.GetPointData().GetScalars().InsertNextValue(length)
self.ExamineSpheres.Modified()
self.vmtkRenderer.RenderWindow.Render()
def Execute(self):
if self.Surface == None:
self.PrintError('Error: No input surface.')
if not self.vmtkRenderer:
self.vmtkRenderer = vmtkrenderer.vmtkRenderer()
self.vmtkRenderer.Initialize()
self.OwnRenderer = 1
self.vmtkRenderer.RegisterScript(self)
glyphs = vtk.vtkGlyph3D()
glyphSource = vtk.vtkSphereSource()
glyphSource.SetRadius(1)
glyphs.SetInput(self.Spheres)
glyphs.SetSource(glyphSource.GetOutput())
glyphs.SetScaleModeToScaleByScalar()
glyphs.SetScaleFactor(1.)
glyphMapper = vtk.vtkPolyDataMapper()
glyphMapper.SetInput(glyphs.GetOutput())
glyphMapper.ScalarVisibilityOff()
self.SpheresActor = vtk.vtkActor()
self.SpheresActor.SetMapper(glyphMapper)
self.SpheresActor.GetProperty().SetColor(1.0,0.0,0.0)
self.SpheresActor.GetProperty().SetOpacity(self.Opacity)
self.SpheresActor.PickableOff()
self.vmtkRenderer.Renderer.AddActor(self.SpheresActor)
examineGlyphs = vtk.vtkGlyph3D()
examineGlyphSource = vtk.vtkSphereSource()
examineGlyphSource.SetRadius(1)
examineGlyphs.SetInput(self.ExamineSpheres)
examineGlyphs.SetSource(examineGlyphSource.GetOutput())
examineGlyphs.SetScaleModeToScaleByScalar()
examineGlyphs.SetScaleFactor(1.)
examineGlyphMapper = vtk.vtkPolyDataMapper()
examineGlyphMapper.SetInput(examineGlyphs.GetOutput())
examineGlyphMapper.ScalarVisibilityOff()
self.ExamineSpheresActor = vtk.vtkActor()
self.ExamineSpheresActor.SetMapper(examineGlyphMapper)
self.ExamineSpheresActor.GetProperty().SetColor(0.0,1.0,0.0)
self.ExamineSpheresActor.GetProperty().SetOpacity(self.Opacity)
self.ExamineSpheresActor.PickableOff()
self.ExamineSpheresActor.VisibilityOff()
self.vmtkRenderer.Renderer.AddActor(self.ExamineSpheresActor)
#self.vmtkRenderer.RenderWindowInteractor.AddObserver("KeyPressEvent", self.KeyPressed)
self.vmtkRenderer.AddKeyBinding('u','undo',self.UndoCallback)
self.vmtkRenderer.AddKeyBinding('plus','Increase sphere radius',self.PlusCallback)
self.vmtkRenderer.AddKeyBinding('minus','Decrease sphere radius',self.MinusCallback)
self.vmtkRenderer.AddKeyBinding('n','Show next sphere',self.NextSphereCallback)
self.vmtkRenderer.AddKeyBinding('v','Show previous sphere',self.PreviousSphereCallback)
self.vmtkRenderer.AddKeyBinding('d','Display ',self.DisplayCallback)
self.vmtkRenderer.AddKeyBinding('w','Switch beetween examien and interact mode ',self.ExmienModeCallback)
self.vmtkRenderer.AddKeyBinding('space','Pick sphere',self.PickCallback)
self.SurfaceMapper = vtk.vtkPolyDataMapper()
self.SurfaceMapper.SetInput(self.Surface)
self.SurfaceMapper.SetScalarVisibility(self.DisplayArray)
surfaceActor = vtk.vtkActor()
surfaceActor.SetMapper(self.SurfaceMapper)
surfaceActor.GetProperty().SetOpacity(self.Opacity)
self.vmtkRenderer.Renderer.AddActor(surfaceActor)
self.ScalarBarActor = vtk.vtkScalarBarActor()
self.ScalarBarActor.SetLookupTable(self.SurfaceMapper.GetLookupTable())
self.ScalarBarActor.GetLabelTextProperty().ItalicOff()
self.ScalarBarActor.GetLabelTextProperty().BoldOff()
self.ScalarBarActor.GetLabelTextProperty().ShadowOff()
self.ScalarBarActor.SetLabelFormat('%.2f')
self.ScalarBarActor.SetTitle('distances')
self.ScalarBarActor.VisibilityOff()
self.vmtkRenderer.Renderer.AddActor(self.ScalarBarActor)
self.SphereWidget = vtk.vtkSphereWidget()
self.SphereWidget.TranslationOff()
self.SphereWidget.SetInteractor(self.vmtkRenderer.RenderWindowInteractor)
self.SphereWidget.AddObserver("InteractionEvent", self.SphereCallback)
self.ExamineText = vtk.vtkTextActor()
self.ExamineText.SetInput("Examine Mode")
self.ExamineText.GetPositionCoordinate().SetCoordinateSystemToNormalizedViewport()
self.ExamineText.SetPosition(0.05,0.95)
self.ExamineText.VisibilityOff()
self.vmtkRenderer.Renderer.AddActor2D(self.ExamineText)
self.InputInfo('Please position the mouse and press space to add spheres, \'u\' to undo\n')
any = 0
while any == 0:
self.InitializeSpheres()
self.vmtkRenderer.Render()
any = (self.Spheres.GetNumberOfPoints()>1)
self.Surface = self.ComputeArray()
if self.Surface.GetSource():
self.Surface.GetSource().UnRegisterAllOutputs()
if self.OwnRenderer:
self.vmtkRenderer.Deallocate()
if __name__=='__main__':
main = pypes.pypeMain()
main.Arguments = sys.argv
main.Execute()
|