-
Notifications
You must be signed in to change notification settings - Fork 0
/
vmtkbranchmetrics.py
120 lines (98 loc) · 5.38 KB
/
vmtkbranchmetrics.py
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
#!/usr/bin/env python
## Program: VMTK
## Module: $RCSfile: vmtkbranchmetrics.py,v $
## Language: Python
## Date: $Date: 2005/09/14 09:48:31 $
## Version: $Revision: 1.4 $
## Copyright (c) Luca Antiga, David Steinman. All rights reserved.
## See LICENSE 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.
from __future__ import absolute_import #NEEDS TO STAY AS TOP LEVEL MODULE FOR Py2-3 COMPATIBILITY
import vtk
from vmtk import vtkvmtk
import sys
from vmtk import pypes
class vmtkBranchMetrics(pypes.pypeScript):
def __init__(self):
pypes.pypeScript.__init__(self)
self.Surface = None
self.Centerlines = None
self.ComputeAbscissaMetric = 1
self.ComputeAngularMetric = 1
self.AbscissasArrayName = 'Abscissas'
self.NormalsArrayName = 'ParallelTransportNormals'
self.GroupIdsArrayName = 'GroupIds'
self.CenterlineIdsArrayName = 'CenterlineIds'
self.TractIdsArrayName = 'TractIds'
self.RadiusArrayName = 'MaximumInscribedSphereRadius'
self.BlankingArrayName = 'Blanking'
self.AngularMetricArrayName = 'AngularMetric'
self.AbscissaMetricArrayName = 'AbscissaMetric'
self.SetScriptName('vmtkbranchmetrics')
self.SetScriptDoc('Takes a centerline and input surface, already split into branches (with centerline attributes calculated) \
and calculates the surface angular and abscissa metric.')
self.SetInputMembers([
['Surface','i','vtkPolyData',1,'','','vmtksurfacereader'],
['Centerlines','centerlines','vtkPolyData',1,'','','vmtksurfacereader'],
['ComputeAbscissaMetric','abscissametric','bool',1],
['ComputeAngularMetric','angularmetric','bool',1],
['AbscissasArrayName','abscissasarray','str',1],
['NormalsArrayName','normalsarray','str',1],
['GroupIdsArrayName','groupidsarray','str',1],
['CenterlineIdsArrayName','centerlineidsarray','str',1],
['TractIdsArrayName','tractidsarray','str',1],
['RadiusArrayName','radiusarray','str',1],
['BlankingArrayName','blankingarray','str',1],
['AngularMetricArrayName','angularmetricarray','str',1],
['AbscissaMetricArrayName','abscissametricarray','str',1]
])
self.SetOutputMembers([
['Surface','o','vtkPolyData',1,'','','vmtksurfacewriter'],
['AngularMetricArrayName','angularmetricarray','str',1],
['AbscissaMetricArrayName','abscissametricarray','str',1]
])
def Execute(self):
if self.Surface == None:
self.PrintError('Error: No input surface.')
if self.Centerlines == None:
self.PrintError('Error: No input centerlines.')
if self.ComputeAngularMetric == 1:
self.PrintLog('Computing angular metric')
angularMetricFilter = vtkvmtk.vtkvmtkPolyDataCenterlineAngularMetricFilter()
angularMetricFilter.SetInputData(self.Surface)
angularMetricFilter.SetMetricArrayName(self.AngularMetricArrayName)
angularMetricFilter.SetGroupIdsArrayName(self.GroupIdsArrayName)
angularMetricFilter.SetCenterlines(self.Centerlines)
angularMetricFilter.SetRadiusArrayName(self.RadiusArrayName)
angularMetricFilter.SetCenterlineNormalsArrayName(self.NormalsArrayName)
angularMetricFilter.SetCenterlineGroupIdsArrayName(self.GroupIdsArrayName)
angularMetricFilter.SetCenterlineTractIdsArrayName(self.TractIdsArrayName)
angularMetricFilter.UseRadiusInformationOff()
angularMetricFilter.IncludeBifurcationsOff()
angularMetricFilter.SetBlankingArrayName(self.BlankingArrayName)
angularMetricFilter.SetCenterlineIdsArrayName(self.CenterlineIdsArrayName)
angularMetricFilter.Update()
self.Surface = angularMetricFilter.GetOutput()
if self.ComputeAbscissaMetric == 1:
self.PrintLog('Computing abscissa metric')
abscissaMetricFilter = vtkvmtk.vtkvmtkPolyDataCenterlineAbscissaMetricFilter()
abscissaMetricFilter.SetInputData(self.Surface)
abscissaMetricFilter.SetMetricArrayName(self.AbscissaMetricArrayName)
abscissaMetricFilter.SetGroupIdsArrayName(self.GroupIdsArrayName)
abscissaMetricFilter.SetCenterlines(self.Centerlines)
abscissaMetricFilter.SetRadiusArrayName(self.RadiusArrayName)
abscissaMetricFilter.SetAbscissasArrayName(self.AbscissasArrayName)
abscissaMetricFilter.SetCenterlineGroupIdsArrayName(self.GroupIdsArrayName)
abscissaMetricFilter.SetCenterlineTractIdsArrayName(self.TractIdsArrayName)
abscissaMetricFilter.UseRadiusInformationOff()
abscissaMetricFilter.IncludeBifurcationsOn()
abscissaMetricFilter.SetBlankingArrayName(self.BlankingArrayName)
abscissaMetricFilter.SetCenterlineIdsArrayName(self.CenterlineIdsArrayName)
abscissaMetricFilter.Update()
self.Surface = abscissaMetricFilter.GetOutput()
if __name__=='__main__':
main = pypes.pypeMain()
main.Arguments = sys.argv
main.Execute()