Skip to content

Contour Tree Alignment

CT bones example Image

Pipeline description

This example first loads a scalar field ensemble from disk as a cinema data base. The ensemble consists of 23 time-dependent scalar fields with 10 time steps each. The ensemble is then filtered for the 23 scalar fields of one fixed time point and read as a vtkMultiBlockDataSet. Here we use the CinemaReader, CinemaQuery and CinemaProductReader filters.

In a pre-processing, the scalar fields are topologically simplified by persistence using the TopologicalSimplificationByPersistence filter. The filter is automatically applied to each member of the MultiBlockDataSet. Then, for each simplified member field, the contour tree is computed using the ContourTree module.

The resulting MultiBlock of contour trees is then used as input for the Contour Tree Alignment filter. This alignment is a super tree of all contour trees and can be seen as a representative of the topology of the whole ensemble. Unfortunately, the vtk object representing the alignment does not have any layout information attached. Therefore, we use the PlanarGraphLayout together with a paraview calculator to compute and apply the layout information.

We now want to check which features of the original scalar fields have been matched onto each other. Therefore, we use the ExtractSeletion filter to extract one vertex and attach its segmentationIDs array to the multi block data set representing the segmentations of the contour trees. We also use a Grid Layout to render the multi block in a comparable fashion (right view, above screenshot).

As a last step, we use the ForEach and EndFor filters to iterate the multi block of segmentations and in each iteration, we extract the region of the scalar field that corresponds to the segmentation id from the selected vertex. The extraction is done using the TTKExtract filter.

ParaView

To reproduce the above screenshot, go to your ttk-data directory and enter the following command:

paraview states/contourTreeAlignment.pvsm

Python code

  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
#!/usr/bin/env python

from paraview.simple import *

# create a new 'TTK CinemaReader'
ttkCinemaReader = TTKCinemaReader(
    DatabasePath="./heatedCylinder/heatedCylinder_2D_raw.cdb"
)

# create a new 'TTK CinemaQuery'
ttkCinemaQuery = TTKCinemaQuery(InputTable=ttkCinemaReader)
ttkCinemaQuery.SQLStatement = """SELECT * FROM InputTable0 WHERE time=4.2"""

# create a new 'TTK CinemaProductReader'
ttkCinemaProductReader = TTKCinemaProductReader(Input=ttkCinemaQuery)

# create a new 'TTK TopologicalSimplificationByPersistence'
ttkTopologicalSimplificationByPersistence = TTKTopologicalSimplificationByPersistence(
    Input=ttkCinemaProductReader
)
ttkTopologicalSimplificationByPersistence.InputArray = ["POINTS", "nrrd"]
ttkTopologicalSimplificationByPersistence.PersistenceThreshold = 0.05
ttkTopologicalSimplificationByPersistence.ThresholdIsAbsolute = True

# create a new 'TTK Merge and Contour Tree ()'
ttkMergeandContourTree = TTKContourTree(
    Input=ttkTopologicalSimplificationByPersistence
)
ttkMergeandContourTree.ScalarField = ["POINTS", "nrrd"]

# create a new 'TTK ContourTreeAlignment'
contourTreeAlignment = TTKContourTreeAlignment(
    Input=OutputPort(ttkMergeandContourTree, 1), ExportPath=""
)
contourTreeAlignment.ScalarField = ["POINTS", "Scalar"]
contourTreeAlignment.Regionsizearray = ["CELLS", "RegionSize"]
contourTreeAlignment.SegmentationIDarrayforCT = ["CELLS", "SegmentationId"]
contourTreeAlignment.SegmentIDarrayforsegmentation = ["POINTS", "Scalar"]
contourTreeAlignment.Seed = 35

# create a new 'TTK PlanarGraphLayout'
alignmentLayout = TTKPlanarGraphLayout(Input=contourTreeAlignment)
alignmentLayout.SequenceArray = ["POINTS", "Scalar"]
alignmentLayout.SizeArray = ["POINTS", "BranchIDs"]
alignmentLayout.UseBranches = 1
alignmentLayout.BranchArray = ["POINTS", "BranchIDs"]
alignmentLayout.LevelArray = ["POINTS", "BranchIDs"]

# create a new 'Calculator'
alignmentEdges = Calculator(Input=alignmentLayout)
alignmentEdges.CoordinateResults = 1
alignmentEdges.Function = "iHat*Layout_Y+jHat*Scalar*3"

# create a query selection
Show()
QuerySelect(
    QueryString="(id == 16)", Source=alignmentEdges, FieldType="POINT", InsideOut=0
)

# create a new 'Extract Selection'
selectedVertex = ExtractSelection(Input=alignmentEdges)

# create a new 'TTK GridLayout'
segmentationsGrid = TTKGridLayout(Input=OutputPort(ttkMergeandContourTree, 2))
segmentationsGrid.ColumnGap = 10.0
segmentationsGrid.RowGap = 10.0

# create a new 'TTK ForEach'
ttkForEach = TTKForEach(Input=segmentationsGrid)
ttkForEach.InputArray = ["POINTS", "SegmentationId"]
ttkForEach.ImageExtent = [0, 127, 0, 255, 0, 0]

# create a new 'Merge Blocks'
mergeBlocks = MergeBlocks(Input=ttkForEach)

# create a new 'TTK ArrayEditor'
passSegmentationIDs = TTKArrayEditor(Target=mergeBlocks, Source=selectedVertex)
passSegmentationIDs.EditorMode = "Add Arrays from Source"
passSegmentationIDs.TargetAttributeType = "Field Data"
passSegmentationIDs.SourcePointDataArrays = ["segmentationIDs"]
passSegmentationIDs.TargetArray = ["POINTS", "SegmentationId"]

# create a new 'TTK Extract'
extractMatchedGeometry = TTKExtract(Input=passSegmentationIDs)
extractMatchedGeometry.ExtractionMode = "Geometry"
extractMatchedGeometry.Expression = "{segmentationIDs[{_ttk_IterationInfo[0]}]}"
extractMatchedGeometry.ImageExtent = [
    2147483647,
    -2147483647,
    2147483647,
    -2147483647,
    2147483647,
    -2147483647,
]
extractMatchedGeometry.InputArray = ["POINTS", "SegmentationId"]

# create a new 'TTK BlockAggregator'
ttkBlockAggregator = TTKBlockAggregator(Input=extractMatchedGeometry)

# create a new 'TTK EndFor'
ttkEndFor = TTKEndFor(Data=ttkBlockAggregator, For=ttkForEach)

# save the output
SaveData("ContourTreeAlignment.vtu", alignmentEdges)
SaveData("Segmentations.vtm", segmentationsGrid)
SaveData("MatchedRegions.vtm", ttkEndFor)

To run the above Python script, go to your ttk-data directory and enter the following command:

pvpython python/contourTreeAlignment.py

Inputs

Outputs

  • ContorTreeAlignment.vtu: the output alignment in VTK file format (left view, above screenshot).
  • Segmentations.vtm: the segmentations of the input scalar fields in VTK multiblock format (right view, above screenshot).
  • MatchedRegions.vtm: the regions of the original fields that are represented by a selected vertex in VTK multiblock format (right view, abÅ¿ove screenshot).

C++/Python API

ArrayEditor

BlockAggregator

CinemaProductReader

CinemaQuery

CinemaReader

ContourTree

ContourTreeAlignment

EndFor

Extract

ForEach

GridLayout

PlanarGraphLayout

TopologicalSimplificationByPersistence