TTK
Loading...
Searching...
No Matches
ttkAlgorithm.cpp
Go to the documentation of this file.
1#include <ttkAlgorithm.h>
2#include <ttkMacros.h>
3#include <ttkUtils.h>
4
6#include <Triangulation.h>
8
9#include <vtkCellTypes.h>
10#include <vtkCommand.h>
11#include <vtkDataSet.h>
12
13#ifdef TTK_ENABLE_MPI
15#include <Identifiers.h>
16#include <vtkCellData.h>
17#include <vtkGhostCellsGenerator.h>
18#endif // TTK_ENABLE_MPI
19
20#include <vtkImageData.h>
21#include <vtkInformation.h>
22#include <vtkInformationIntegerKey.h>
23#include <vtkInformationVector.h>
24#include <vtkMultiBlockDataSet.h>
25#include <vtkPointData.h>
26#include <vtkPolyData.h>
27#include <vtkTable.h>
28#include <vtkUnstructuredGrid.h>
29
30#include <vtkCompositeDataPipeline.h>
31
32// Pass input type information key
33#include <vtkInformationKey.h>
34vtkInformationKeyMacro(ttkAlgorithm, SAME_DATA_TYPE_AS_INPUT_PORT, Integer);
35
36// Constructor / Destructor
40
42
43 this->printMsg("Requesting triangulation for '"
44 + std::string(dataSet->GetClassName()) + "'",
46#ifdef TTK_ENABLE_MPI
47 if((ttk::hasInitializedMPI()) && (ttk::isRunningWithMPI())) {
48 if(!hasMPISupport_) {
49 printErr("MPI is not formally supported for this filter :(");
50 printErr("The results are likely to be incorrect.");
51 }
53 }
54#endif // TTK_ENABLE_MPI
55
57 this->debugLevel_, this->CompactTriangulationCacheSize, dataSet);
58
59#ifdef TTK_ENABLE_MPI
60 if(ttk::hasInitializedMPI()) {
61 std::vector<int> tmp{};
62 std::map<int, int> tmpId{};
63 this->MPIPipelinePreconditioning(dataSet, tmp, tmpId, triangulation);
64 this->MPITriangulationPreconditioning(triangulation, dataSet);
65 }
66#endif // TTK_ENABLE_MPI
67
68 if(triangulation)
69 return triangulation;
70
71 this->printErr("Unable to retrieve/initialize triangulation for '"
72 + std::string(dataSet->GetClassName()) + "'");
73
74 return nullptr;
75}
76
77vtkDataArray *ttkAlgorithm::GetOptionalArray(const bool &enforceArrayIndex,
78 const int &arrayIndex,
79 const std::string &arrayName,
80 vtkDataSet *const inputData,
81 const int &inputPort) {
82
83 vtkDataArray *optionalArray = nullptr;
84
85 if(enforceArrayIndex)
86 optionalArray = this->GetInputArrayToProcess(arrayIndex, inputData);
87
88 if(!optionalArray) {
89 this->SetInputArrayToProcess(arrayIndex, inputPort, 0, 0, arrayName.data());
90 optionalArray = this->GetInputArrayToProcess(arrayIndex, inputData);
91 }
92 return optionalArray;
93}
94
95std::string ttkAlgorithm::GetOrderArrayName(vtkDataArray *const array) {
96 return std::string(array->GetName()) + "_Order";
97}
98
99vtkDataArray *
100 ttkAlgorithm::ComputeOrderArray(vtkDataSet *const inputData,
101 vtkDataArray *scalarArray,
102 const int scalarArrayIdx,
103 const bool getGlobalOrder,
104 vtkDataArray *oldOrderArray,
105 ttk::Triangulation *triangulation) {
106
108 auto nVertices = scalarArray->GetNumberOfTuples();
109 if(oldOrderArray != nullptr && getGlobalOrder) {
110 newOrderArray = ttkSimplexIdTypeArray::SafeDownCast(oldOrderArray);
111 } else {
113 newOrderArray->SetName(this->GetOrderArrayName(scalarArray).data());
114 newOrderArray->SetNumberOfComponents(1);
115 newOrderArray->SetNumberOfTuples(nVertices);
116 }
117
118 std::vector<int> neighbors;
119 std::map<int, int> neighborsToId;
120#ifdef TTK_ENABLE_MPI
121 if(ttk::hasInitializedMPI()) {
122 this->MPIGhostPipelinePreconditioning(inputData);
124 inputData, neighbors, neighborsToId, nullptr);
125 }
126 if(ttk::isRunningWithMPI() && getGlobalOrder) {
127 ttk::ArrayPreconditioning arrayPreconditioning
129 arrayPreconditioning.preconditionTriangulation(triangulation);
130 arrayPreconditioning.setGlobalOrder(getGlobalOrder);
131 ttkTypeMacroAT(scalarArray->GetDataType(), triangulation->getType(),
132 (arrayPreconditioning.processScalarArray<T0, T1>(
133 static_cast<const T1 *>(triangulation->getData()),
135 ttkUtils::GetPointer<T0>(scalarArray), nVertices)));
136 } else {
137 switch(scalarArray->GetDataType()) {
138 vtkTemplateMacro(ttk::preconditionOrderArray(
139 nVertices, static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(scalarArray)),
140 static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(newOrderArray)),
141 this->threadNumber_));
142 }
143 }
144 if(oldOrderArray == nullptr || !getGlobalOrder) {
145 inputData
146 ->GetAttributesAsFieldData(
147 this->GetInputArrayAssociation(scalarArrayIdx, inputData))
148 ->AddArray(newOrderArray);
149 }
150#else
151
152 switch(scalarArray->GetDataType()) {
153 vtkTemplateMacro(ttk::preconditionOrderArray(
154 nVertices, static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(scalarArray)),
155 static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(newOrderArray)),
156 this->threadNumber_));
157 }
158
159 inputData
160 ->GetAttributesAsFieldData(
161 this->GetInputArrayAssociation(scalarArrayIdx, inputData))
162 ->AddArray(newOrderArray);
163 TTK_FORCE_USE(triangulation);
164#endif
165 return newOrderArray;
166}
167
169 vtkDataSet *const inputData,
170 vtkDataArray *scalarArray,
171 const int scalarArrayIdx,
172 const bool getGlobalOrder,
173 vtkDataArray *orderArray,
174 ttk::Triangulation *triangulation,
175 const bool enforceOrderArrayIdx) {
176
177 std::string enforcedArray = "";
178 if(enforceOrderArrayIdx) {
179 enforcedArray = " enforced ";
180 }
181#ifdef TTK_ENABLE_MPI
182 if(getGlobalOrder) {
183 if(triangulation->isOrderArrayGlobal(
184 ttkUtils::GetVoidPointer(scalarArray))) {
185 this->printMsg("Retrieved " + enforcedArray + " order array `"
186 + std::string(orderArray->GetName()) + "`.",
188 return orderArray;
189 } else {
190 ttk::Timer timer;
192 this->printWrn("Order array `" + std::string(orderArray->GetName())
193 + "` is local, but a global order array is "
194 "required. Re-computing.");
195
196 this->printMsg("Initializing order array.", 0, 0, this->threadNumber_,
199
200 orderArray
201 = this->ComputeOrderArray(inputData, scalarArray, scalarArrayIdx,
202 getGlobalOrder, orderArray, triangulation);
203
204 triangulation->setIsOrderArrayGlobal(
205 ttkUtils::GetVoidPointer(scalarArray), true);
206
207 this->printMsg("Initializing order array.", 1, timer.getElapsedTime(),
208 this->threadNumber_);
209
211 this->printWrn("TIP: run `ttkArrayPreconditioning` first with "
212 "GlobalOrder enabled");
213 this->printWrn("for improved performances :)");
215 return orderArray;
216 }
217 } else {
218#else
219 TTK_FORCE_USE(inputData);
220 TTK_FORCE_USE(scalarArray);
221 TTK_FORCE_USE(scalarArrayIdx);
222 TTK_FORCE_USE(getGlobalOrder);
223 TTK_FORCE_USE(triangulation);
224#endif // TTK_ENABLE_MPI
225 this->printMsg("Retrieved " + enforcedArray + " order array `"
226 + std::string(orderArray->GetName()) + "`.",
228 return orderArray;
229#ifdef TTK_ENABLE_MPI
230 }
231#endif
232}
233
234vtkDataArray *ttkAlgorithm::GetOrderArray(vtkDataSet *const inputData,
235 const int scalarArrayIdx,
236 ttk::Triangulation *triangulation,
237 const bool getGlobalOrder,
238 const int orderArrayIdx,
239 const bool enforceOrderArrayIdx) {
240
241 auto isValidOrderArray = [](vtkDataArray *const array) {
242 if(!array)
243 return -4;
244
245 if(array->GetNumberOfComponents() != 1)
246 return -3;
247
249 if(array->GetDataType() != temp->GetDataType())
250 return -2;
251
252 const std::string name(array->GetName());
253 if(name.size() < 6 || (name.rfind("_Order") != (name.size() - 6)))
254 return -1;
255
256 return 1;
257 };
258 auto scalarArray = this->GetInputArrayToProcess(scalarArrayIdx, inputData);
259 if(enforceOrderArrayIdx) {
260 auto orderArray = this->GetInputArrayToProcess(orderArrayIdx, inputData);
261 switch(isValidOrderArray(orderArray)) {
262 case -4: {
263 this->printErr("Unable to retrieve enforced order array at idx "
264 + std::to_string(orderArrayIdx) + ".");
265 return nullptr;
266 }
267 case -3: {
268 this->printErr("Retrieved enforced order array `"
269 + std::string(orderArray->GetName())
270 + "` has more than one component.");
271 return nullptr;
272 }
273 case -2: {
274 this->printErr("Enforced order array `"
275 + std::string(orderArray->GetName())
276 + "` is of incorrect type.");
278 this->printErr(" -> use `ttkArrayEditor` to convert data type to `"
279 + std::string(temp->GetDataTypeAsString()) + "`.");
280 return nullptr;
281 }
282 default: {
284 inputData, scalarArray, scalarArrayIdx, getGlobalOrder, orderArray,
285 triangulation, enforceOrderArrayIdx);
286 }
287 }
288 }
289
290 if(!scalarArray) {
291 this->printErr("Unable to retrieve input scalar array for idx "
292 + std::to_string(scalarArrayIdx) + ".");
293 return nullptr;
294 } else if(isValidOrderArray(scalarArray) == 1) {
295 this->printMsg("Retrieved scalar array `"
296 + std::string(scalarArray->GetName())
297 + "` is already an order array.",
299 return scalarArray;
300 }
301
302 auto orderArray = inputData
303 ->GetAttributesAsFieldData(this->GetInputArrayAssociation(
304 scalarArrayIdx, inputData))
305 ->GetArray(this->GetOrderArrayName(scalarArray).data());
306
307 switch(isValidOrderArray(orderArray)) {
308 case -4: {
309 ttk::Timer timer;
311 this->printWrn("No pre-existing order for array:");
312 this->printWrn(" `" + std::string(scalarArray->GetName()) + "`.");
313
314 this->printMsg("Initializing order array.", 0, 0, this->threadNumber_,
317
318 orderArray
319 = this->ComputeOrderArray(inputData, scalarArray, scalarArrayIdx,
320 getGlobalOrder, orderArray, triangulation);
321
322 std::string optionOn = "";
323#ifdef TTK_ENABLE_MPI
324 if(getGlobalOrder) {
325 optionOn = "with GlobalOrder enabled ";
326 }
327 bool isGlobalOrder = getGlobalOrder || (!ttk::isRunningWithMPI());
328 triangulation->setIsOrderArrayGlobal(
329 ttkUtils::GetVoidPointer(scalarArray), isGlobalOrder);
330#endif // TTK_ENABLE_MPI
331 this->printMsg("Initializing order array.", 1, timer.getElapsedTime(),
332 this->threadNumber_);
333
335 this->printWrn("TIP: run `ttkArrayPreconditioning` first");
336 this->printWrn(optionOn + "for improved performances :)");
338
339 return orderArray;
340 }
341
342 case -3: {
343 this->printErr(
344 "Retrieved order array `" + std::string(orderArray->GetName())
345 + "` for scalar array `" + std::string(scalarArray->GetName())
346 + "` has more than one component.");
347 return nullptr;
348 }
349
350 case -2: {
351 this->printErr(
352 "Retrieved order array `" + std::string(orderArray->GetName())
353 + "` for scalar array `" + std::string(scalarArray->GetName())
354 + "` is of incorrect type.");
356 this->printErr(" -> use `ttkArrayEditor` to convert data type to `"
357 + std::string(temp->GetDataTypeAsString()) + "`.");
358 return nullptr;
359 }
360
361 default: {
362#ifdef TTK_ENABLE_MPI
363 if(!ttk::isRunningWithMPI()) {
364#endif
365 this->printMsg("Retrieved order array `"
366 + std::string(orderArray->GetName())
367 + "` for scalar array `"
368 + std::string(scalarArray->GetName()) + "`.",
370 return orderArray;
371#ifdef TTK_ENABLE_MPI
372 }
374 inputData, scalarArray, scalarArrayIdx, getGlobalOrder, orderArray,
375 triangulation, enforceOrderArrayIdx);
376#endif
377 }
378 }
379}
380
382 ttkAlgorithm::GetIdentifierArrayPtr(const bool &enforceArrayIndex,
383 const int &arrayIndex,
384 const std::string &arrayName,
385 vtkDataSet *const inputData,
386 std::vector<ttk::SimplexId> &spareStorage,
387 const int inputPort,
388 const bool printErr) {
389
390 // fetch data array
391 const auto array = this->GetOptionalArray(
392 enforceArrayIndex, arrayIndex, arrayName, inputData, inputPort);
393 if(array == nullptr) {
394 if(printErr) {
395 this->printErr("Could not find the requested identifiers array");
396 }
397 return {};
398 }
399 if(array->GetNumberOfComponents() != 1) {
400 if(printErr) {
401 this->printErr("Identifiers field must have only one component!");
402 }
403 return {};
404 }
405
406#ifndef TTK_ENABLE_64BIT_IDS
407 if(array->GetDataType() == VTK_ID_TYPE
408 || array->GetDataType() == VTK_LONG_LONG) {
409 this->printMsg(
410 "Converting identifiers field from vtkIdType to SimplexId...");
411 const auto nItems = array->GetNumberOfTuples();
412
413 // fills the vector with the content of the data array converted to
414 // ttk::SimplexId
415 spareStorage.resize(nItems);
416 for(vtkIdType i = 0; i < nItems; ++i) {
417 spareStorage[i] = static_cast<ttk::SimplexId>(array->GetTuple1(i));
418 }
419
420 // return a pointer to the vector internal buffer
421 return spareStorage.data();
422 }
423#else
424 TTK_FORCE_USE(spareStorage);
425#endif // TTK_ENABLE_64BIT_IDS
426
427 // return a pointer to the data array internal buffer
428 return static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(array));
429}
430
431template <class vtkDataType>
432int prepOutput(vtkInformation *info, const std::string &className) {
433 auto output = vtkDataObject::GetData(info);
434 if(!output || !output->IsA(className.data())) {
435 auto newOutput = vtkSmartPointer<vtkDataType>::New();
436 info->Set(vtkDataObject::DATA_OBJECT(), newOutput);
437 }
438 return 1;
439}
440
442 return this->GetOutput(0);
443}
444
445vtkDataSet *ttkAlgorithm::GetOutput(int port) {
446 return vtkDataSet::SafeDownCast(this->GetOutputDataObject(port));
447}
448
449void ttkAlgorithm::SetInputData(vtkDataSet *input) {
450 this->SetInputData(0, input);
451}
452
453void ttkAlgorithm::SetInputData(int index, vtkDataSet *input) {
454 this->SetInputDataInternal(index, input);
455}
456
457void ttkAlgorithm::AddInputData(vtkDataSet *input) {
458 this->AddInputData(0, input);
459}
460
461void ttkAlgorithm::AddInputData(int index, vtkDataSet *input) {
462 this->AddInputDataInternal(index, input);
463}
464
465int ttkAlgorithm::RequestDataObject(vtkInformation *ttkNotUsed(request),
466 vtkInformationVector **inputVector,
467 vtkInformationVector *outputVector) {
468 // for each output
469 for(int i = 0; i < this->GetNumberOfOutputPorts(); ++i) {
470 auto outInfo = outputVector->GetInformationObject(i);
471 if(!outInfo) {
472 this->printErr("Unable to retrieve output vtkDataObject at port "
473 + std::to_string(i));
474 return 0;
475 }
476
477 auto outputPortInfo = this->GetOutputPortInformation(i);
478
479 // always request output type again for dynamic filter outputs
480 if(!this->FillOutputPortInformation(i, outputPortInfo)) {
481 this->printErr("Unable to fill output port information at port "
482 + std::to_string(i));
483 return 0;
484 }
485
486 if(outputPortInfo->Has(ttkAlgorithm::SAME_DATA_TYPE_AS_INPUT_PORT())) {
487 // Set output data type to input data type at specified port
488 auto inPortIndex
489 = outputPortInfo->Get(ttkAlgorithm::SAME_DATA_TYPE_AS_INPUT_PORT());
490 if(inPortIndex < 0 || inPortIndex >= this->GetNumberOfInputPorts()) {
491 this->printErr("Input port index " + std::to_string(inPortIndex)
492 + " specified by 'SAME_DATA_TYPE_AS_INPUT_PORT' key of "
493 "output port is out of range ("
494 + std::to_string(this->GetNumberOfInputPorts())
495 + " input ports).");
496 return 0;
497 }
498 auto inInfo = inputVector[inPortIndex]->GetInformationObject(0);
499 if(!inInfo) {
500 this->printErr(
501 "No information object at port " + std::to_string(inPortIndex)
502 + " specified by 'SAME_DATA_TYPE_AS_INPUT_PORT' key of output port.");
503 return 0;
504 }
505
506 auto input = vtkDataObject::GetData(inInfo);
507 auto output = vtkDataObject::GetData(outInfo);
508
509 if(!output || !output->IsA(input->GetClassName())) {
510 auto newOutput
511 = vtkSmartPointer<vtkDataObject>::Take(input->NewInstance());
512 outputPortInfo->Set(
513 vtkDataObject::DATA_TYPE_NAME(), input->GetClassName());
514 outInfo->Set(vtkDataObject::DATA_OBJECT(), newOutput);
515 }
516 } else {
517 // Explicitly create output by data type name
518 if(!outputPortInfo->Has(vtkDataObject::DATA_TYPE_NAME())) {
519 this->printErr("DATA_TYPE_NAME of output port " + std::to_string(i)
520 + " not specified");
521 return 0;
522 }
523 std::string const outputType
524 = outputPortInfo->Get(vtkDataObject::DATA_TYPE_NAME());
525
526 if(outputType == "vtkUnstructuredGrid") {
527 prepOutput<vtkUnstructuredGrid>(outInfo, outputType);
528 } else if(outputType == "vtkPolyData") {
529 prepOutput<vtkPolyData>(outInfo, outputType);
530 } else if(outputType == "vtkMultiBlockDataSet") {
531 prepOutput<vtkMultiBlockDataSet>(outInfo, outputType);
532 } else if(outputType == "vtkTable") {
533 prepOutput<vtkTable>(outInfo, outputType);
534 } else if(outputType == "vtkImageData") {
535 prepOutput<vtkImageData>(outInfo, outputType);
536 } else {
537 this->printErr("Unsupported data type for output[" + std::to_string(i)
538 + "]: " + outputType);
539 return 0;
540 }
541 }
542
543 this->printMsg(
544 "Created '"
545 + std::string(outputPortInfo->Get(vtkDataObject::DATA_TYPE_NAME()))
546 + "' at output port " + std::to_string(i),
548 }
549
550 return 1;
551}
552
553#ifdef TTK_ENABLE_MPI
554
555int ttkAlgorithm::updateMPICommunicator(vtkDataSet *input) {
556 if(input == nullptr) {
557 return 0;
558 }
559 int isEmpty
560 = input->GetNumberOfCells() == 0 || input->GetNumberOfPoints() == 0;
561 int oldSize = ttk::MPIsize_;
562 int oldRank = ttk::MPIrank_;
563 MPI_Comm_split(MPI_COMM_WORLD, isEmpty, 0, &ttk::MPIcomm_);
564 MPI_Comm_rank(ttk::MPIcomm_, &ttk::MPIrank_);
565 MPI_Comm_size(ttk::MPIcomm_, &ttk::MPIsize_);
566 if(oldSize != ttk::MPIsize_) {
567 std::vector<int> newToOldRanks(ttk::MPIsize_);
568 MPI_Allgather(&oldRank, 1, MPI_INTEGER, newToOldRanks.data(), 1,
569 MPI_INTEGER, ttk::MPIcomm_);
570 std::map<int, int> oldToNewRanks;
571 for(int i = 0; i < ttk::MPIsize_; i++) {
572 oldToNewRanks[newToOldRanks[i]] = i;
573 }
574 int *vertexRankArray
575 = ttkUtils::GetPointer<int>(input->GetPointData()->GetArray("RankArray"));
576 if(vertexRankArray != nullptr) {
577 for(int i = 0; i < input->GetNumberOfPoints(); i++) {
578 vertexRankArray[i] = oldToNewRanks[vertexRankArray[i]];
579 }
580 }
581 int *cellRankArray
582 = ttkUtils::GetPointer<int>(input->GetCellData()->GetArray("RankArray"));
583 if(cellRankArray != nullptr) {
584 for(int i = 0; i < input->GetNumberOfCells(); i++) {
585 cellRankArray[i] = oldToNewRanks[cellRankArray[i]];
586 }
587 }
588 }
590 return isEmpty;
591}
592
594 ttk::SimplexId simplexNumber,
595 unsigned char *ghost,
596 int *rankArray) {
597 ttk::SimplexId ghostNumber = 0;
598 if(rankArray != nullptr) {
599#ifdef TTK_ENABLE_OPENMP
600#pragma omp parallel for reduction(+ : ghostNumber)
601#endif // TTK_ENABLE_OPENMP
602 for(ttk::SimplexId i = 0; i < simplexNumber; i++) {
603 if(rankArray[i] != ttk::MPIrank_) {
604 ghostNumber++;
605 }
606 }
607 } else {
608 if(ghost != nullptr) {
609#ifdef TTK_ENABLE_OPENMP
610#pragma omp parallel for reduction(+ : ghostNumber)
611#endif // TTK_ENABLE_OPENMP
612 for(ttk::SimplexId i = 0; i < simplexNumber; i++) {
613 if(ghost[i] == 1) {
614 ghostNumber++;
615 }
616 }
617 }
618 }
619
620 ttk::SimplexId realSimplexNumber = simplexNumber - ghostNumber;
621 auto minmax = std::minmax_element(globalIds, globalIds + simplexNumber);
622 ttk::LongSimplexId min = globalIds[minmax.first - globalIds];
623 ttk::LongSimplexId max = globalIds[minmax.second - globalIds];
624 ttk::SimplexId globalSimplexNumber;
625 ttk::LongSimplexId globalMin;
626 ttk::LongSimplexId globalMax;
627 MPI_Allreduce(&realSimplexNumber, &globalSimplexNumber, 1,
628 ttk::getMPIType(realSimplexNumber), MPI_SUM, ttk::MPIcomm_);
629 MPI_Allreduce(
630 &min, &globalMin, 1, ttk::getMPIType(min), MPI_MIN, ttk::MPIcomm_);
631 MPI_Allreduce(
632 &max, &globalMax, 1, ttk::getMPIType(max), MPI_MAX, ttk::MPIcomm_);
633
634 return (globalSimplexNumber == globalMax + 1 && globalMin == 0);
635};
636
638 vtkDataSet *input,
639 std::unordered_map<ttk::SimplexId, ttk::SimplexId> &vertGtoL,
640 std::vector<int> &neighborRanks,
641 std::map<int, int> &neighborsToId) {
642
643 ttk::Identifiers identifiers;
644
645 vtkNew<vtkIdTypeArray> vtkVertexIdentifiers{};
646
647 vtkNew<vtkIdTypeArray> vtkCellIdentifiers{};
648 vtkVertexIdentifiers->SetName("GlobalPointIds");
649 vtkVertexIdentifiers->SetNumberOfComponents(1);
650 vtkVertexIdentifiers->SetNumberOfTuples(input->GetNumberOfPoints());
651 vtkVertexIdentifiers->Fill(-1);
652
653 identifiers.setVertexIdentifiers(
654 ttkUtils::GetPointer<ttk::LongSimplexId>(vtkVertexIdentifiers));
655
656 vtkCellIdentifiers->SetName("GlobalCellIds");
657 vtkCellIdentifiers->SetNumberOfComponents(1);
658 vtkCellIdentifiers->SetNumberOfTuples(input->GetNumberOfCells());
659 vtkCellIdentifiers->Fill(-1);
660
661 identifiers.setCellIdentifiers(
662 ttkUtils::GetPointer<ttk::LongSimplexId>(vtkCellIdentifiers));
663
664 int vertexNumber = input->GetNumberOfPoints();
665 identifiers.setVertexNumber(vertexNumber);
666 int cellNumber = input->GetNumberOfCells();
667 identifiers.setCellNumber(cellNumber);
668 int status = 0;
669
670 double *boundingBox = input->GetBounds();
671 identifiers.setBounds(boundingBox);
672 identifiers.initializeNeighbors(boundingBox, neighborRanks, neighborsToId);
673 if(ttk::isRunningWithMPI()) {
674 switch(input->GetDataObjectType()) {
675 case VTK_UNSTRUCTURED_GRID:
676 case VTK_POLY_DATA: {
677
678 identifiers.setOutdatedGlobalPointIds(
680 input->GetPointData()->GetGlobalIds()));
681 identifiers.setOutdatedGlobalCellIds(
683 input->GetCellData()->GetGlobalIds()));
684 identifiers.setVertexRankArray(ttkUtils::GetPointer<int>(
685 input->GetPointData()->GetArray("RankArray")));
686 int *cellRankArray = ttkUtils::GetPointer<int>(
687 input->GetCellData()->GetArray("RankArray"));
688 identifiers.setCellRankArray(cellRankArray);
689 identifiers.setVertGhost(ttkUtils::GetPointer<unsigned char>(
690 input->GetPointData()->GetArray("vtkGhostType")));
691 unsigned char *cellGhost = ttkUtils::GetPointer<unsigned char>(
692 input->GetCellData()->GetArray("vtkGhostType"));
693 identifiers.setCellGhost(cellGhost);
694 vtkPointSet *pointSet = vtkPointSet::SafeDownCast(input);
695 identifiers.setPointSet(static_cast<float *>(
696 ttkUtils::GetVoidPointer(pointSet->GetPoints())));
697 vtkCellArray *cells = nullptr;
698 switch(input->GetDataObjectType()) {
699 case VTK_UNSTRUCTURED_GRID: {
700 auto dataSetAsUG = vtkUnstructuredGrid::SafeDownCast(input);
701 cells = dataSetAsUG->GetCells();
702 break;
703 }
704 case VTK_POLY_DATA: {
705 auto dataSetAsPD = vtkPolyData::SafeDownCast(input);
706 cells
707 = dataSetAsPD->GetNumberOfPolys() > 0 ? dataSetAsPD->GetPolys()
708 : dataSetAsPD->GetNumberOfLines() > 0 ? dataSetAsPD->GetLines()
709 : dataSetAsPD->GetVerts();
710 break;
711 }
712 default: {
713 this->printErr("Unable to get cells for `"
714 + std::string(input->GetClassName()) + "`");
715 }
716 }
717 if(cells == nullptr) {
718 return 0;
719 }
720 if(!cells->IsStorage64Bit()) {
721 if(cells->CanConvertTo64BitStorage()) {
722 this->printWrn("Converting the cell array to 64-bit storage");
723 bool success = cells->ConvertTo64BitStorage();
724 if(!success) {
725 this->printErr(
726 "Error converting the provided cell array to 64-bit storage");
727 return -1;
728 }
729 } else {
730 this->printErr(
731 "Cannot convert the provided cell array to 64-bit storage");
732 return -1;
733 }
734 }
735
736 identifiers.setConnectivity(ttkUtils::GetPointer<ttk::LongSimplexId>(
737 cells->GetConnectivityArray()));
738
739 std::vector<std::vector<ttk::SimplexId>> pointsToCells(vertexNumber);
740 vtkIdList *cellList = vtkIdList::New();
741 if(cellRankArray != nullptr) {
742 for(ttk::SimplexId i = 0; i < vertexNumber; i++) {
743 input->GetPointCells(i, cellList);
744 for(int j = 0; j < cellList->GetNumberOfIds(); j++) {
745 if(cellRankArray[cellList->GetId(j)] == ttk::MPIrank_) {
746 pointsToCells[i].push_back(cellList->GetId(j));
747 }
748 }
749 }
750 } else {
751 for(ttk::SimplexId i = 0; i < vertexNumber; i++) {
752 input->GetPointCells(i, cellList);
753 for(int j = 0; j < cellList->GetNumberOfIds(); j++) {
754 if(cellGhost[cellList->GetId(j)] == 0) {
755 pointsToCells[i].push_back(cellList->GetId(j));
756 }
757 }
758 }
759 }
760 identifiers.setPointsToCells(pointsToCells);
761
762 identifiers.initializeMPITypes();
763 identifiers.setVertGtoL(&vertGtoL);
764 vtkIdList *pointCell = vtkIdList::New();
765 input->GetCellPoints(0, pointCell);
766 int nbPoints = pointCell->GetNumberOfIds();
767 identifiers.setDomainDimension(nbPoints - 1);
768 identifiers.buildKDTree();
769 status = identifiers.executePolyData();
770 break;
771 }
772 case VTK_IMAGE_DATA: {
773 vtkImageData *data = vtkImageData::SafeDownCast(input);
774 identifiers.setDims(data->GetDimensions());
775 identifiers.setSpacing(data->GetSpacing());
776 status = identifiers.executeImageData();
777 break;
778 }
779 default: {
780 this->printErr("Unable to triangulate `"
781 + std::string(input->GetClassName()) + "`");
782 }
783 }
784 } else {
785 status = identifiers.executeSequential();
786 }
787
788 if(status < 1) {
789 printErr("Global identifier generation failed");
790 return -1;
791 }
792
793 // Add VTK objects to the data set
794
795 input->GetPointData()->SetGlobalIds(vtkVertexIdentifiers);
796 input->GetCellData()->SetGlobalIds(vtkCellIdentifiers);
797 return 0;
798}
799
800void ttkAlgorithm::MPIGhostPipelinePreconditioning(vtkDataSet *input) {
801
802 vtkNew<vtkGhostCellsGenerator> generator;
803 if(ttk::isRunningWithMPI()
804 && (!input->HasAnyGhostCells()
805 && ((input->GetPointData()->GetArray("RankArray") == nullptr)
806 || (input->GetCellData()->GetArray("RankArray") == nullptr)))) {
807 generator->SetInputData(input);
808 generator->BuildIfRequiredOff();
809 generator->SetNumberOfGhostLayers(1);
810 generator->Update();
811 input->ShallowCopy(generator->GetOutputDataObject(0));
812 input->GetPointData()->AddArray(
813 generator->GetOutputDataObject(0)->GetGhostArray(0));
814 input->GetCellData()->AddArray(
815 generator->GetOutputDataObject(0)->GetGhostArray(1));
816 }
817}
818
820 vtkDataSet *input,
821 std::vector<int> &neighbors,
822 std::map<int, int> &neighToId,
823 ttk::Triangulation *triangulation) {
824
825 ttk::SimplexId vertexNumber = input->GetNumberOfPoints();
826 ttk::SimplexId cellNumber = input->GetNumberOfCells();
827
828 if((input->GetDataObjectType() == VTK_POLY_DATA
829 || input->GetDataObjectType() == VTK_UNSTRUCTURED_GRID)) {
830
831 if((ttk::hasInitializedMPI()) && (ttk::isRunningWithMPI())) {
833 printWrn("The distribution by VTK of Unstructured");
834 printWrn("Grids and Poly Data has been reported");
835 printWrn("to be affected by bugs (at least up");
836 printWrn("to ParaView 5.10.1).");
837
838 if((input->GetCellData()->GetGlobalIds() == nullptr)
839 || (input->GetPointData()->GetGlobalIds() == nullptr)) {
840
841 printWrn("=> Global identifiers may be incorrect.");
842 }
843 if((input->GetPointData()->GetArray("RankArray") == nullptr)
844 || (input->GetCellData()->GetArray("RankArray") == nullptr)) {
845 printWrn("=> Rank arrays may be incorrect.");
846 }
848 }
849 }
850
851 // Get the neighbor ranks
852 std::vector<int> &neighborRanks{
853 triangulation != nullptr ? triangulation->getNeighborRanks() : neighbors};
854 std::map<int, int> &neighborsToId{
855 triangulation != nullptr ? triangulation->getNeighborsToId() : neighToId};
856
857 double *boundingBox = input->GetBounds();
858 if(triangulation != nullptr) {
859 triangulation->createMetaGrid(boundingBox);
860 }
861
862 if(neighborRanks.empty()) {
863 ttk::preconditionNeighborsUsingBoundingBox(
864 boundingBox, neighborRanks, neighborsToId);
865 }
866
867 // Checks if global ids are valid
869 input->GetPointData()->GetGlobalIds());
871 input->GetCellData()->GetGlobalIds());
872
873 bool pointValidity{false};
874 bool cellValidity{false};
875 if((triangulation != nullptr
876 && (triangulation->getType() == ttk::Triangulation::Type::EXPLICIT
877 || triangulation->getType() == ttk::Triangulation::Type::COMPACT))
878 || triangulation == nullptr) {
879 if(globalPointIds != nullptr) {
880 unsigned char *ghostPoints = ttkUtils::GetPointer<unsigned char>(
881 input->GetPointData()->GetArray("vtkGhostType"));
882 int *vertexRankArray = ttkUtils::GetPointer<int>(
883 input->GetPointData()->GetArray("RankArray"));
884 pointValidity = checkGlobalIdValidity(
885 globalPointIds, vertexNumber, ghostPoints, vertexRankArray);
886 }
887 if(pointValidity && globalCellIds != nullptr) {
888
889 unsigned char *ghostCells = ttkUtils::GetPointer<unsigned char>(
890 input->GetCellData()->GetArray("vtkGhostType"));
891 int *cellRankArray = ttkUtils::GetPointer<int>(
892 input->GetCellData()->GetArray("RankArray"));
893 cellValidity = checkGlobalIdValidity(
894 globalCellIds, cellNumber, ghostCells, cellRankArray);
895 }
896 } else {
897 pointValidity = true;
898 cellValidity = true;
899 }
900
901 // If the global ids are not valid, they are computed again
902 if(!pointValidity || !cellValidity) {
903 if(triangulation != nullptr) {
904 if(triangulation->getType() == ttk::Triangulation::Type::EXPLICIT
905 || triangulation->getType() == ttk::Triangulation::Type::COMPACT) {
906 this->GenerateGlobalIds(input, triangulation->getVertexGlobalIdMap(),
907 neighborRanks, neighborsToId);
908 }
909 } else {
910 std::unordered_map<ttk::SimplexId, ttk::SimplexId> vertGtoL{};
911 this->GenerateGlobalIds(input, vertGtoL, neighborRanks, neighborsToId);
912 }
913 }
914}
915
917 ttk::Triangulation *triangulation, vtkDataSet *input) {
918
919 const auto pd{input->GetPointData()};
920 if(pd == nullptr) {
921 triangulation->printWrn("No point data on input object");
922 } else {
923 // provide "GlobalPointIds" & "vtkGhostType" point data arrays
924 // to the triangulation
925 triangulation->setVertsGlobalIds(
926 ttkUtils::GetPointer<ttk::LongSimplexId>(pd->GetGlobalIds()));
927 triangulation->setVertexGhostArray(
928 ttkUtils::GetPointer<unsigned char>(pd->GetArray("vtkGhostType")));
929 int *vertexRankArray = ttkUtils::GetPointer<int>(pd->GetArray("RankArray"));
930 if(vertexRankArray != nullptr) {
931 triangulation->setVertexRankArray(vertexRankArray);
932 }
933 triangulation->preconditionDistributedVertices();
934 }
935
936 const auto cd{input->GetCellData()};
937 if(cd == nullptr) {
938 triangulation->printWrn("No cell data on input object");
939 } else {
940 // provide "GlobalCellIds" & "vtkGhostType" cell data arrays to
941 // the triangulation
942 triangulation->setCellsGlobalIds(
943 ttkUtils::GetPointer<ttk::LongSimplexId>(cd->GetGlobalIds()));
944 triangulation->setCellGhostArray(
945 ttkUtils::GetPointer<unsigned char>(cd->GetArray("vtkGhostType")));
946 int *cellRankArray = ttkUtils::GetPointer<int>(cd->GetArray("RankArray"));
947 if(cellRankArray != nullptr) {
948 triangulation->setCellRankArray(cellRankArray);
949 }
950 triangulation->preconditionDistributedCells();
951 }
952}
953
954#endif // TTK_ENABLE_MPI
955
956//==============================================================================
957int ttkAlgorithm::ProcessRequest(vtkInformation *request,
958 vtkInformationVector **inputVector,
959 vtkInformationVector *outputVector) {
960 // 1. Pass
961 if(request->Has(vtkCompositeDataPipeline::REQUEST_DATA_OBJECT())) {
962 this->printMsg(
963 "Processing REQUEST_DATA_OBJECT", ttk::debug::Priority::VERBOSE);
964 return this->RequestDataObject(request, inputVector, outputVector);
965 }
966
967 // 2. Pass
968 if(request->Has(vtkCompositeDataPipeline::REQUEST_INFORMATION())) {
969 this->printMsg(
970 "Processing REQUEST_INFORMATION", ttk::debug::Priority::VERBOSE);
971 return this->RequestInformation(request, inputVector, outputVector);
972 }
973
974 // 3. Pass
975 if(request->Has(vtkCompositeDataPipeline::REQUEST_UPDATE_TIME())) {
976 this->printMsg(
977 "Processing REQUEST_UPDATE_TIME", ttk::debug::Priority::VERBOSE);
978 return this->RequestUpdateTime(request, inputVector, outputVector);
979 }
980
981 // 4. Pass
982 if(request->Has(
983 vtkCompositeDataPipeline::REQUEST_TIME_DEPENDENT_INFORMATION())) {
984 this->printMsg("Processing REQUEST_TIME_DEPENDENT_INFORMATION",
987 request, inputVector, outputVector);
988 }
989
990 // 5. Pass
991 if(request->Has(vtkCompositeDataPipeline::REQUEST_UPDATE_EXTENT())) {
992 this->printMsg(
993 "Processing REQUEST_UPDATE_EXTENT", ttk::debug::Priority::VERBOSE);
994 return this->RequestUpdateExtent(request, inputVector, outputVector);
995 }
996
997 // 6. Pass
998 if(request->Has(vtkCompositeDataPipeline::REQUEST_DATA_NOT_GENERATED())) {
999 this->printMsg(
1000 "Processing REQUEST_DATA_NOT_GENERATED", ttk::debug::Priority::VERBOSE);
1001 return this->RequestDataNotGenerated(request, inputVector, outputVector);
1002 }
1003
1004 // 7. Pass
1005 if(request->Has(vtkCompositeDataPipeline::REQUEST_DATA())) {
1006 this->printMsg("Processing REQUEST_DATA", ttk::debug::Priority::VERBOSE);
1008#ifdef TTK_ENABLE_MPI
1009 if(ttk::hasInitializedMPI() && inputVector != nullptr) {
1010 if(this->updateMPICommunicator(vtkDataSet::GetData(inputVector[0], 0))) {
1011 return 1;
1012 };
1013 }
1014#endif // TTK_ENABLE_MPI
1015 return this->RequestData(request, inputVector, outputVector);
1016 }
1017
1018 this->printErr("Unsupported pipeline pass:");
1019 request->Print(cout);
1020
1021 return 0;
1022}
#define TTK_FORCE_USE(x)
Force the compiler to use the function/method parameter.
Definition BaseClass.h:57
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
Baseclass of all VTK filters that wrap ttk modules.
void SetInputData(vtkDataSet *)
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
bool checkGlobalIdValidity(ttk::LongSimplexId *globalIds, ttk::SimplexId simplexNumber, unsigned char *ghost, int *rankArray)
ttk::SimplexId * GetIdentifierArrayPtr(const bool &enforceArrayIndex, const int &arrayIndex, const std::string &arrayName, vtkDataSet *const inputData, std::vector< ttk::SimplexId > &spareStorage, const int inputPort=0, const bool printErr=true)
void AddInputData(vtkDataSet *)
virtual int RequestUpdateTime(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
~ttkAlgorithm() override
virtual int RequestDataNotGenerated(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
float CompactTriangulationCacheSize
virtual int RequestData(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
void MPIGhostPipelinePreconditioning(vtkDataSet *input)
virtual int RequestInformation(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
int ProcessRequest(vtkInformation *request, vtkInformationVector **inputVectors, vtkInformationVector *outputVector) override
int GenerateGlobalIds(vtkDataSet *input, std::unordered_map< ttk::SimplexId, ttk::SimplexId > &vertGtoL, std::vector< int > &neighborRanks, std::map< int, int > &neighborsToId)
vtkDataArray * checkForGlobalAndComputeOrderArray(vtkDataSet *const inputData, vtkDataArray *scalarArray, const int scalarArrayIdx, const bool getGlobalOrder, vtkDataArray *orderArray, ttk::Triangulation *triangulation, const bool enforceOrderArrayIdx)
virtual int RequestDataObject(vtkInformation *request, vtkInformationVector **inputVectors, vtkInformationVector *outputVector)
virtual int RequestUpdateTimeDependentInformation(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
virtual int RequestUpdateExtent(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
void MPITriangulationPreconditioning(ttk::Triangulation *triangulation, vtkDataSet *input)
vtkDataSet * GetOutput()
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, ttk::Triangulation *triangulation, const bool getGlobalOrder=false, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
void MPIPipelinePreconditioning(vtkDataSet *input, std::vector< int > &neighbors, std::map< int, int > &neighToId, ttk::Triangulation *triangulation=nullptr)
int FillOutputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
vtkDataArray * ComputeOrderArray(vtkDataSet *const inputData, vtkDataArray *scalarArray, const int scalarArrayIdx, const bool getGlobalOrder, vtkDataArray *oldOrderArray, ttk::Triangulation *triangulation)
vtkDataArray * GetOptionalArray(const bool &enforceArrayIndex, const int &arrayIndex, const std::string &arrayName, vtkDataSet *const inputData, const int &inputPort=0)
static std::string GetOrderArrayName(vtkDataArray *const array)
static ttk::Triangulation * GetTriangulation(int debugLevel, float cacheRatio, vtkDataSet *object)
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
static DT * GetPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.h:59
int preconditionTriangulation(AbstractTriangulation *triangulation)
int processScalarArray(const triangulationType *triangulation, ttk::SimplexId *orderArray, const DT *scalarArray, const size_t nVerts) const
int debugLevel_
Definition Debug.h:379
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
std::string debugMsgNamePrefix_
Definition Debug.h:384
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void setCellIdentifiers(ttk::LongSimplexId *cellIdentifiers)
Definition Identifiers.h:97
void setVertexNumber(const SimplexId &vertexNumber)
Definition Identifiers.h:85
void setCellNumber(const SimplexId &cellNumber)
Definition Identifiers.h:89
int executeSequential()
Generates global ids for all data set type in sequential.
void setVertexIdentifiers(ttk::LongSimplexId *vertexIdentifiers)
Definition Identifiers.h:93
double getElapsedTime()
Definition Timer.h:15
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
Triangulation::Type getType() const
FiltratedEdge max(const FiltratedEdge &a, const FiltratedEdge &b)
COMMON_EXPORTS int MPIsize_
Definition BaseClass.cpp:10
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
COMMON_EXPORTS int MPIrank_
Definition BaseClass.cpp:9
void preconditionOrderArray(const size_t nVerts, const scalarType *const scalars, SimplexId *const order, const int nThreads=ttk::globalThreadNumber_)
Precondition an order array to be consumed by the base layer API.
long long int LongSimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:15
vtkStandardNewMacro(ttkAlgorithm)
vtkInformationKeyMacro(ttkAlgorithm, SAME_DATA_TYPE_AS_INPUT_PORT, Integer)
int prepOutput(vtkInformation *info, const std::string &className)
#define ttkTypeMacroAT(group0, group1, call)
Definition ttkMacros.h:272
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)