TTK
Loading...
Searching...
No Matches
ttkConnectedComponents.cpp
Go to the documentation of this file.
2
3#include <vtkDataArray.h>
4#include <vtkInformation.h>
5#include <vtkObjectFactory.h>
6#include <vtkPointData.h>
7#include <vtkPolyData.h>
8#include <vtkSmartPointer.h>
9
10#include <vtkFloatArray.h>
11#include <vtkIntArray.h>
12
13#include <ttkMacros.h>
14#include <ttkUtils.h>
15
17
19 this->SetNumberOfInputPorts(1);
20 this->SetNumberOfOutputPorts(2);
21
22 // Suppress warning if one does not set the optional input array
23 this->SetInputArrayToProcess(0, 0, 0, 0, "%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%");
24}
25
27 vtkInformation *info) {
28 if(port == 0) {
29 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
30 return 1;
31 }
32 return 0;
33}
34
36 vtkInformation *info) {
37 if(port == 0) {
39 return 1;
40 } else if(port == 1) {
41 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData");
42 return 1;
43 }
44 return 0;
45}
46
48 vtkInformationVector **inputVector,
49 vtkInformationVector *outputVector) {
50
51 // Fetch Input Data
52 auto inputDataSet = vtkDataSet::GetData(inputVector[0]);
53 if(!inputDataSet)
54 return 0;
55
56 const auto nPoints = inputDataSet->GetNumberOfPoints();
57
58 auto featureMask = this->GetInputArrayToProcess(0, inputVector);
59
60 if(featureMask && this->GetInputArrayAssociation(0, inputVector) != 0)
61 return !this->printErr("Input array needs to be a point data array.");
62
63 if(featureMask && featureMask->GetNumberOfComponents() != 1)
64 return !this->printErr("Input array needs to be a scalar array.");
65
66 auto triangulation = ttkAlgorithm::GetTriangulation(inputDataSet);
67 if(!triangulation)
68 return 0;
69
70 // Allocate Output Label Data
71 auto componentIds = vtkSmartPointer<vtkIntArray>::New();
72 componentIds->SetName("ComponentId");
73 componentIds->SetNumberOfTuples(nPoints);
74
75 auto componentSize = vtkSmartPointer<vtkFloatArray>::New();
76 componentSize->SetName("ComponentSize");
77 if(this->AugmentSegmentationWithComponentSize)
78 componentSize->SetNumberOfTuples(nPoints);
79
80 // Compute Connected Components
81 std::vector<ttk::ConnectedComponents::Component> components;
82 {
83 int status = 0;
84
85 // initialize component ids
86 ttkTypeMacroA(featureMask ? featureMask->GetDataType() : VTK_INT,
87 (status = this->initializeComponentIds<T0>(
88 ttkUtils::GetPointer<int>(componentIds), nPoints,
89 ttkUtils::GetPointer<const T0>(featureMask),
90 this->BackgroundThreshold)));
91 if(!status)
92 return 0;
93
94 // compute connected components (componentIds now store component idx)
95 this->preconditionTriangulation(triangulation);
96 ttkTypeMacroT(triangulation->getType(),
97 (status = this->computeConnectedComponents<T0>(
98 components, ttkUtils::GetPointer<int>(componentIds),
99 static_cast<const T0 *>(triangulation->getData()))));
100 if(!status)
101 return 0;
102 }
103
104 // Augment Segmentation
105 {
106 int status = 0;
107
108 // map component size to segmentation
109 if(this->AugmentSegmentationWithComponentSize) {
110 status = this->mapData(
111 components, ttkUtils::GetPointer<int>(componentIds), nPoints,
112 [](const std::vector<ttk::ConnectedComponents::Component> &components_,
113 const int cIdx) { return cIdx < 0 ? 0.0 : components_[cIdx].size; },
114 ttkUtils::GetPointer<float>(componentSize));
115 if(!status)
116 return 0;
117 }
118
119 // map component idx to component id
120 status = this->mapData(
121 components, ttkUtils::GetPointer<int>(componentIds), nPoints,
122 [](const std::vector<ttk::ConnectedComponents::Component> &components_,
123 const int cIdx) { return cIdx < 0 ? -1 : components_[cIdx].id; },
124 ttkUtils::GetPointer<int>(componentIds));
125 if(!status)
126 return 0;
127 }
128
129 // Segmentation Output
130 {
131 auto outputDataSet = vtkDataSet::GetData(outputVector, 0);
132 outputDataSet->ShallowCopy(inputDataSet);
133 outputDataSet->GetPointData()->AddArray(componentIds);
134 if(this->AugmentSegmentationWithComponentSize)
135 outputDataSet->GetPointData()->AddArray(componentSize);
136 }
137
138 // Components Output
139 {
140 const int nComponents = components.size();
141 auto outputComponents = vtkPolyData::GetData(outputVector, 1);
142
143 // points
144 {
145 auto sizeArray = vtkSmartPointer<vtkFloatArray>::New();
146 sizeArray->SetName("Size");
147 sizeArray->SetNumberOfTuples(nComponents);
148 auto sizeArrayData = ttkUtils::GetPointer<float>(sizeArray);
149
150 auto idArray = vtkSmartPointer<vtkIntArray>::New();
151 idArray->SetName("ComponentId");
152 idArray->SetNumberOfTuples(nComponents);
153 auto idArrayData = ttkUtils::GetPointer<int>(idArray);
154
155 auto points = vtkSmartPointer<vtkPoints>::New();
156 points->SetDataTypeToFloat();
157 points->SetNumberOfPoints(nComponents);
158 auto pointsData = ttkUtils::GetPointer<float>(points->GetData());
159 for(int i = 0, j = 0; i < nComponents; i++) {
160 const auto &c = components[i];
161 pointsData[j++] = c.center[0];
162 pointsData[j++] = c.center[1];
163 pointsData[j++] = c.center[2];
164
165 sizeArrayData[i] = c.size;
166 idArrayData[i] = c.id;
167 }
168
169 outputComponents->SetPoints(points);
170 auto pd = outputComponents->GetPointData();
171 pd->AddArray(sizeArray);
172 pd->AddArray(idArray);
173 }
174
175 // cells
176 {
177 auto connectivityArray = vtkSmartPointer<vtkIntArray>::New();
178 connectivityArray->SetNumberOfTuples(nComponents);
179 auto connectivityArrayData = ttkUtils::GetPointer<int>(connectivityArray);
180 for(int i = 0; i < nComponents; i++)
181 connectivityArrayData[i] = i;
182
183 auto offsetArray = vtkSmartPointer<vtkIntArray>::New();
184 offsetArray->SetNumberOfTuples(nComponents + 1);
185 auto offsetArrayData = ttkUtils::GetPointer<int>(offsetArray);
186 for(int i = 0; i <= nComponents; i++)
187 offsetArrayData[i] = i;
188
189 auto cellArray = vtkSmartPointer<vtkCellArray>::New();
190 cellArray->SetData(offsetArray, connectivityArray);
191
192 outputComponents->SetVerts(cellArray);
193 }
194
195 // Copy Field Data
196 outputComponents->GetFieldData()->ShallowCopy(inputDataSet->GetFieldData());
197 }
198
199 // return success
200 return 1;
201}
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter that computes connected components based on a scalar field.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
int preconditionTriangulation(ttk::AbstractTriangulation *triangulation) const
int mapData(std::vector< Component > &components, const int *componentIds, const int nVertices, const F f, DT *out) const
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
vtkStandardNewMacro(ttkConnectedComponents)
#define ttkTypeMacroT(group, call)
Definition ttkMacros.h:195
#define ttkTypeMacroA(group, call)
Definition ttkMacros.h:234