TTK
Loading...
Searching...
No Matches
ttkMergeTreeDistanceMatrix.cpp
Go to the documentation of this file.
2#include <ttkMergeTreeUtils.h>
3#include <ttkUtils.h>
4
5#include <vtkDataObject.h> // For port information
6#include <vtkObjectFactory.h> // for new macro
7
8#include <vtkDoubleArray.h>
9#include <vtkInformation.h>
10#include <vtkStringArray.h>
11#include <vtkTable.h>
12
13using namespace ttk;
14using namespace ftm;
15
16// A VTK macro that enables the instantiation of this class via ::New()
17// You do not have to modify this
19
33 this->SetNumberOfInputPorts(2);
34 this->SetNumberOfOutputPorts(1);
35}
36
38
47 vtkInformation *info) {
48 if(port == 0 || port == 1) {
49 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
50 if(port == 1)
51 info->Set(vtkAlgorithm::INPUT_IS_OPTIONAL(), 1);
52 } else
53 return 0;
54
55 return 1;
56}
57
74 int port, vtkInformation *info) {
75 if(port == 0)
76 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkTable");
77 else
78 return 0;
79
80 return 1;
81}
82
96template <class dataType>
98 vtkInformationVector *outputVector,
99 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
100 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2) {
101
102 // Construct trees
103 const int numInputs = inputTrees.size();
104 std::vector<MergeTree<dataType>> intermediateTrees, intermediateTrees2;
105 bool const useSecondPairsType
106 = (mixtureCoefficient_ == 0); // only for PD support
108 inputTrees, intermediateTrees, useSecondPairsType, DiagramPairTypes);
110 or (mixtureCoefficient_ != 0 and mixtureCoefficient_ != 1)) {
111 auto &inputTrees2ToUse
112 = (not isPersistenceDiagram_ ? inputTrees2 : inputTrees);
113 constructTrees(inputTrees2ToUse, intermediateTrees2, !useSecondPairsType,
114 DiagramPairTypes);
115 }
116
117 // Verify parameters
118 if(not UseFieldDataParameters) {
119 if(Backend == 0) {
122 keepSubtree_ = false;
123 baseModule_ = 0;
124 } else if(Backend == 1) {
125 branchDecomposition_ = false;
127 keepSubtree_ = true;
128 baseModule_ = 0;
129 } else if(Backend == 3) {
132 keepSubtree_ = true;
133 baseModule_ = 1;
134 } else if(Backend == 4) {
137 keepSubtree_ = true;
138 baseModule_ = 2;
139 } else {
140 baseModule_ = 0;
141 }
142 }
143 if(baseModule_ == 0) {
146 }
147 if(not branchDecomposition_) {
149 printMsg("NormalizedWasserstein is set to false since branch "
150 "decomposition is not asked.");
152 }
154 printMsg("Computation with normalized Wasserstein.");
155 else
156 printMsg("Computation without normalized Wasserstein.");
160 printMsg("BranchDecomposition: " + std::to_string(branchDecomposition_));
161 printMsg("NormalizedWasserstein: "
162 + std::to_string(normalizedWasserstein_));
163 printMsg("KeepSubtree: " + std::to_string(keepSubtree_));
164 }
165 if(baseModule_ == 1) {
166 printMsg("Using Branch Mapping Distance.");
167 std::string metric;
168 if(branchMetric_ == 0)
169 metric = "Wasserstein Distance first degree";
170 else if(branchMetric_ == 1)
171 metric = "Wasserstein Distance second degree";
172 else if(branchMetric_ == 2)
173 metric = "Persistence difference";
174 else if(branchMetric_ == 3)
175 metric = "Shifting cost";
176 else
177 return 1;
178 printMsg("BranchMetric: " + metric);
179 }
180 if(baseModule_ == 2) {
181 printMsg("Using Path Mapping Distance.");
182 std::string metric;
183 if(pathMetric_ == 0)
184 metric = "Persistence difference";
185 else
186 return 1;
187 printMsg("PathMetric: " + metric);
188 }
189
190 // --- Call base
191 std::vector<std::vector<double>> treesDistMat(
192 numInputs, std::vector<double>(numInputs));
193 if(baseModule_ == 0)
194 execute<dataType>(intermediateTrees, intermediateTrees2, treesDistMat);
195 else
196 execute<dataType>(intermediateTrees, treesDistMat);
197
198 // --- Create output
199 auto treesDistTable = vtkTable::GetData(outputVector);
200
201 // zero-padd column name to keep Row Data columns ordered
202 const auto zeroPad
203 = [](std::string &colName, const size_t numberCols, const size_t colIdx) {
204 std::string const max{std::to_string(numberCols - 1)};
205 std::string const cur{std::to_string(colIdx)};
206 std::string const zer(max.size() - cur.size(), '0');
207 colName.append(zer).append(cur);
208 };
209
210 // copy trees distance matrix to output
211 vtkNew<vtkIntArray> treeIds{};
212 treeIds->SetName("treeID");
213 treeIds->SetNumberOfTuples(numInputs);
214 for(size_t i = 0; i < treesDistMat.size(); ++i) {
215 treeIds->SetTuple1(i, i);
216
217 std::string name{"Tree"};
218 zeroPad(name, treesDistMat.size(), i);
219 vtkNew<vtkDoubleArray> col{};
220 col->SetNumberOfTuples(numInputs);
221 col->SetName(name.c_str());
222 for(size_t j = 0; j < treesDistMat[i].size(); ++j) {
223 col->SetTuple1(j, treesDistMat[i][j]);
224 }
225 treesDistTable->AddColumn(col);
226 }
227
228 treesDistTable->AddColumn(treeIds);
229
230 // aggregate input field data
231 vtkNew<vtkFieldData> allFieldData{}, allFieldDataCopy{};
232 for(unsigned int i = 0; i < inputTrees.size(); ++i) {
233 for(unsigned int j = 0; j < inputTrees[i]->GetNumberOfBlocks(); ++j) {
234 auto fd = inputTrees[i]->GetBlock(j)->GetFieldData();
235 for(int k = 0; k < fd->GetNumberOfArrays(); ++k) {
236 auto array = fd->GetAbstractArray(k);
237 auto dataArray = vtkDataArray::SafeDownCast(array);
238 auto stringArray = vtkStringArray::SafeDownCast(array);
239 if(dataArray or stringArray)
240 allFieldData->AddArray(array);
241 }
242 }
243 }
244 allFieldDataCopy->DeepCopy(allFieldData); // to not modify original field data
245
246 for(int k = 0; k < allFieldDataCopy->GetNumberOfArrays(); ++k) {
247 auto array = allFieldDataCopy->GetAbstractArray(k);
248 array->SetNumberOfTuples(inputTrees.size());
249 auto dataArray = vtkDataArray::SafeDownCast(array);
250 auto stringArray = vtkStringArray::SafeDownCast(array);
251 auto name = array->GetName();
252 for(unsigned int i = 0; i < inputTrees.size(); ++i) {
253 bool foundArray = false;
254 for(unsigned int j = 0; j < inputTrees[i]->GetNumberOfBlocks(); ++j) {
255 auto inputArray
256 = inputTrees[i]->GetBlock(j)->GetFieldData()->GetAbstractArray(name);
257 if(inputArray) {
258 array->SetTuple(i, 0, inputArray);
259 foundArray = true;
260 } else if(not foundArray) {
261 if(dataArray) {
262 const double val = std::nan("");
263 dataArray->SetTuple(i, &val);
264 } else if(stringArray) {
265 stringArray->SetValue(i, "");
266 }
267 }
268 }
269 }
270 treesDistTable->AddColumn(array);
271 }
272
273 return 1;
274}
275
277 vtkInformation *ttkNotUsed(request),
278 vtkInformationVector **inputVector,
279 vtkInformationVector *outputVector) {
280 // --- Get input object from input vector
281 auto blocks = vtkMultiBlockDataSet::GetData(inputVector[0], 0);
282 auto blocks2 = vtkMultiBlockDataSet::GetData(inputVector[1], 0);
283
284 // --- Load blocks
285 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> inputTrees, inputTrees2;
286 loadBlocks(inputTrees, blocks);
287 loadBlocks(inputTrees2, blocks2);
288
289 // --- Load field data parameters
290 if(UseFieldDataParameters) {
291 printMsg("Load parameters from field data.");
292 std::vector<std::string> paramNames;
293 getParamNames(paramNames);
294 for(auto paramName : paramNames) {
295 auto array = blocks->GetFieldData()->GetArray(paramName.c_str());
296 if(array) {
297 double const value = array->GetTuple1(0);
298 setParamValueFromName(paramName, value);
299 printMsg(" - " + paramName + " = " + std::to_string(value));
300 } else
301 printMsg(" - " + paramName + " was not found in the field data.");
302 }
303 }
304
305 return run<float>(outputVector, inputTrees, inputTrees2);
306}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
TTK VTK-filter that wraps the ttk::MergeTreeDistanceMatrix module.
int FillInputPortInformation(int port, vtkInformation *info) override
~ttkMergeTreeDistanceMatrix() override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int run(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
void setParamValueFromName(std::string &paramName, double value)
void getParamNames(std::vector< std::string > &paramNames)
void execute(std::vector< ftm::MergeTree< dataType > > &trees, std::vector< ftm::MergeTree< dataType > > &trees2, std::vector< std::vector< double > > &distanceMatrix)
bool constructTrees(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< MergeTree< dataType > > &intermediateTrees, std::vector< vtkUnstructuredGrid * > &treesNodes, std::vector< vtkUnstructuredGrid * > &treesArcs, std::vector< vtkDataSet * > &treesSegmentation, const std::vector< bool > &useSecondPairsTypeVec, int diagramPairTypes=0)
void loadBlocks(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, vtkMultiBlockDataSet *blocks)
TTK base package defining the standard types.
vtkStandardNewMacro(ttkMergeTreeDistanceMatrix)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)