9#include <vtkCellTypes.h>
10#include <vtkCommand.h>
11#include <vtkDataSet.h>
16#include <vtkCellData.h>
17#include <vtkGhostCellsGenerator.h>
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>
28#include <vtkUnstructuredGrid.h>
30#include <vtkCompositeDataPipeline.h>
33#include <vtkInformationKey.h>
43 this->
printMsg(
"Requesting triangulation for '"
44 + std::string(dataSet->GetClassName()) +
"'",
47 if((ttk::hasInitializedMPI()) && (ttk::isRunningWithMPI())) {
49 printErr(
"MPI is not formally supported for this filter :(");
50 printErr(
"The results are likely to be incorrect.");
60 if(ttk::hasInitializedMPI()) {
61 std::vector<int> tmp{};
62 std::map<int, int> tmpId{};
71 this->
printErr(
"Unable to retrieve/initialize triangulation for '"
72 + std::string(dataSet->GetClassName()) +
"'");
78 const int &arrayIndex,
79 const std::string &arrayName,
80 vtkDataSet *
const inputData,
81 const int &inputPort) {
83 vtkDataArray *optionalArray =
nullptr;
86 optionalArray = this->GetInputArrayToProcess(arrayIndex, inputData);
89 this->SetInputArrayToProcess(arrayIndex, inputPort, 0, 0, arrayName.data());
90 optionalArray = this->GetInputArrayToProcess(arrayIndex, inputData);
96 return std::string(array->GetName()) +
"_Order";
101 vtkDataArray *scalarArray,
102 const int scalarArrayIdx,
103 const bool getGlobalOrder,
104 vtkDataArray *oldOrderArray,
108 auto nVertices = scalarArray->GetNumberOfTuples();
109 if(oldOrderArray !=
nullptr && getGlobalOrder) {
110 newOrderArray = ttkSimplexIdTypeArray::SafeDownCast(oldOrderArray);
114 newOrderArray->SetNumberOfComponents(1);
115 newOrderArray->SetNumberOfTuples(nVertices);
118 std::vector<int> neighbors;
119 std::map<int, int> neighborsToId;
121 if(ttk::hasInitializedMPI()) {
124 inputData, neighbors, neighborsToId,
nullptr);
126 if(ttk::isRunningWithMPI() && getGlobalOrder) {
133 static_cast<const T1 *
>(triangulation->
getData()),
134 ttkUtils::GetPointer<ttk::SimplexId>(newOrderArray),
135 ttkUtils::GetPointer<T0>(scalarArray), nVertices)));
137 switch(scalarArray->GetDataType()) {
144 if(oldOrderArray ==
nullptr || !getGlobalOrder) {
146 ->GetAttributesAsFieldData(
147 this->GetInputArrayAssociation(scalarArrayIdx, inputData))
148 ->AddArray(newOrderArray);
151 switch(scalarArray->GetDataType()) {
158 ->GetAttributesAsFieldData(
159 this->GetInputArrayAssociation(scalarArrayIdx, inputData))
160 ->AddArray(newOrderArray);
163 return newOrderArray;
167 vtkDataSet *
const inputData,
168 vtkDataArray *scalarArray,
169 const int scalarArrayIdx,
170 const bool getGlobalOrder,
171 vtkDataArray *orderArray,
173 const bool enforceOrderArrayIdx) {
175 std::string enforcedArray =
"";
176 if(enforceOrderArrayIdx) {
177 enforcedArray =
" enforced ";
181 if(triangulation->isOrderArrayGlobal(
183 this->
printMsg(
"Retrieved " + enforcedArray +
" order array `"
184 + std::string(orderArray->GetName()) +
"`.",
190 this->
printWrn(
"Order array `" + std::string(orderArray->GetName())
191 +
"` is local, but a global order array is "
192 "required. Re-computing.");
200 getGlobalOrder, orderArray, triangulation);
202 triangulation->setIsOrderArrayGlobal(
206 this->threadNumber_);
209 this->
printWrn(
"TIP: run `ttkArrayPreconditioning` first with "
210 "GlobalOrder enabled");
211 this->
printWrn(
"for improved performances :)");
223 this->
printMsg(
"Retrieved " + enforcedArray +
" order array `"
224 + std::string(orderArray->GetName()) +
"`.",
233 const int scalarArrayIdx,
235 const bool getGlobalOrder,
236 const int orderArrayIdx,
237 const bool enforceOrderArrayIdx) {
239 auto isValidOrderArray = [](vtkDataArray *
const array) {
243 if(array->GetNumberOfComponents() != 1)
247 if(array->GetDataType() != temp->GetDataType())
250 const std::string name(array->GetName());
251 if(name.size() < 6 || (name.rfind(
"_Order") != (name.size() - 6)))
256 auto scalarArray = this->GetInputArrayToProcess(scalarArrayIdx, inputData);
257 if(enforceOrderArrayIdx) {
258 auto orderArray = this->GetInputArrayToProcess(orderArrayIdx, inputData);
259 switch(isValidOrderArray(orderArray)) {
261 this->
printErr(
"Unable to retrieve enforced order array at idx "
262 + std::to_string(orderArrayIdx) +
".");
266 this->
printErr(
"Retrieved enforced order array `"
267 + std::string(orderArray->GetName())
268 +
"` has more than one component.");
272 this->
printErr(
"Enforced order array `"
273 + std::string(orderArray->GetName())
274 +
"` is of incorrect type.");
276 this->
printErr(
" -> use `ttkArrayEditor` to convert data type to `"
277 + std::string(temp->GetDataTypeAsString()) +
"`.");
282 inputData, scalarArray, scalarArrayIdx, getGlobalOrder, orderArray,
283 triangulation, enforceOrderArrayIdx);
289 this->
printErr(
"Unable to retrieve input scalar array for idx "
290 + std::to_string(scalarArrayIdx) +
".");
292 }
else if(isValidOrderArray(scalarArray) == 1) {
293 this->
printMsg(
"Retrieved scalar array `"
294 + std::string(scalarArray->GetName())
295 +
"` is already an order array.",
300 auto orderArray = inputData
301 ->GetAttributesAsFieldData(this->GetInputArrayAssociation(
302 scalarArrayIdx, inputData))
305 switch(isValidOrderArray(orderArray)) {
309 this->
printWrn(
"No pre-existing order for array:");
310 this->
printWrn(
" `" + std::string(scalarArray->GetName()) +
"`.");
318 getGlobalOrder, orderArray, triangulation);
320 std::string optionOn =
"";
323 optionOn =
"with GlobalOrder enabled ";
325 bool isGlobalOrder = getGlobalOrder || (!ttk::isRunningWithMPI());
326 triangulation->setIsOrderArrayGlobal(
330 this->threadNumber_);
333 this->
printWrn(
"TIP: run `ttkArrayPreconditioning` first");
334 this->
printWrn(optionOn +
"for improved performances :)");
342 "Retrieved order array `" + std::string(orderArray->GetName())
343 +
"` for scalar array `" + std::string(scalarArray->GetName())
344 +
"` has more than one component.");
350 "Retrieved order array `" + std::string(orderArray->GetName())
351 +
"` for scalar array `" + std::string(scalarArray->GetName())
352 +
"` is of incorrect type.");
354 this->
printErr(
" -> use `ttkArrayEditor` to convert data type to `"
355 + std::string(temp->GetDataTypeAsString()) +
"`.");
361 inputData, scalarArray, scalarArrayIdx, getGlobalOrder, orderArray,
362 triangulation, enforceOrderArrayIdx);
369 const int &arrayIndex,
370 const std::string &arrayName,
371 vtkDataSet *
const inputData,
372 std::vector<ttk::SimplexId> &spareStorage,
374 const bool printErr) {
378 enforceArrayIndex, arrayIndex, arrayName, inputData, inputPort);
379 if(array ==
nullptr) {
381 this->
printErr(
"Could not find the requested identifiers array");
385 if(array->GetNumberOfComponents() != 1) {
387 this->
printErr(
"Identifiers field must have only one component!");
392#ifndef TTK_ENABLE_64BIT_IDS
393 if(array->GetDataType() == VTK_ID_TYPE
394 || array->GetDataType() == VTK_LONG_LONG) {
396 "Converting identifiers field from vtkIdType to SimplexId...");
397 const auto nItems = array->GetNumberOfTuples();
401 spareStorage.resize(nItems);
402 for(vtkIdType i = 0; i < nItems; ++i) {
403 spareStorage[i] =
static_cast<ttk::SimplexId>(array->GetTuple1(i));
407 return spareStorage.data();
417template <
class vtkDataType>
418int prepOutput(vtkInformation *info,
const std::string &className) {
419 auto output = vtkDataObject::GetData(info);
420 if(!output || !output->IsA(className.data())) {
422 info->Set(vtkDataObject::DATA_OBJECT(), newOutput);
432 return vtkDataSet::SafeDownCast(this->GetOutputDataObject(port));
440 this->SetInputDataInternal(index, input);
448 this->AddInputDataInternal(index, input);
452 vtkInformationVector **inputVector,
453 vtkInformationVector *outputVector) {
455 for(
int i = 0; i < this->GetNumberOfOutputPorts(); ++i) {
456 auto outInfo = outputVector->GetInformationObject(i);
458 this->
printErr(
"Unable to retrieve output vtkDataObject at port "
459 + std::to_string(i));
463 auto outputPortInfo = this->GetOutputPortInformation(i);
467 this->
printErr(
"Unable to fill output port information at port "
468 + std::to_string(i));
476 if(inPortIndex < 0 || inPortIndex >= this->GetNumberOfInputPorts()) {
477 this->
printErr(
"Input port index " + std::to_string(inPortIndex)
478 +
" specified by 'SAME_DATA_TYPE_AS_INPUT_PORT' key of "
479 "output port is out of range ("
480 + std::to_string(this->GetNumberOfInputPorts())
484 auto inInfo = inputVector[inPortIndex]->GetInformationObject(0);
487 "No information object at port " + std::to_string(inPortIndex)
488 +
" specified by 'SAME_DATA_TYPE_AS_INPUT_PORT' key of output port.");
492 auto input = vtkDataObject::GetData(inInfo);
493 auto output = vtkDataObject::GetData(outInfo);
495 if(!output || !output->IsA(input->GetClassName())) {
499 vtkDataObject::DATA_TYPE_NAME(), input->GetClassName());
500 outInfo->Set(vtkDataObject::DATA_OBJECT(), newOutput);
504 if(!outputPortInfo->Has(vtkDataObject::DATA_TYPE_NAME())) {
505 this->
printErr(
"DATA_TYPE_NAME of output port " + std::to_string(i)
509 std::string
const outputType
510 = outputPortInfo->Get(vtkDataObject::DATA_TYPE_NAME());
512 if(outputType ==
"vtkUnstructuredGrid") {
513 prepOutput<vtkUnstructuredGrid>(outInfo, outputType);
514 }
else if(outputType ==
"vtkPolyData") {
515 prepOutput<vtkPolyData>(outInfo, outputType);
516 }
else if(outputType ==
"vtkMultiBlockDataSet") {
517 prepOutput<vtkMultiBlockDataSet>(outInfo, outputType);
518 }
else if(outputType ==
"vtkTable") {
519 prepOutput<vtkTable>(outInfo, outputType);
520 }
else if(outputType ==
"vtkImageData") {
521 prepOutput<vtkImageData>(outInfo, outputType);
523 this->
printErr(
"Unsupported data type for output[" + std::to_string(i)
524 +
"]: " + outputType);
531 + std::string(outputPortInfo->Get(vtkDataObject::DATA_TYPE_NAME()))
532 +
"' at output port " + std::to_string(i),
541int ttkAlgorithm::updateMPICommunicator(vtkDataSet *input) {
542 if(input ==
nullptr) {
546 = input->GetNumberOfCells() == 0 || input->GetNumberOfPoints() == 0;
549 MPI_Comm_split(MPI_COMM_WORLD, isEmpty, 0, &ttk::MPIcomm_);
554 MPI_Allgather(&oldRank, 1, MPI_INTEGER, newToOldRanks.data(), 1,
555 MPI_INTEGER, ttk::MPIcomm_);
556 std::map<int, int> oldToNewRanks;
558 oldToNewRanks[newToOldRanks[i]] = i;
561 = ttkUtils::GetPointer<int>(input->GetPointData()->GetArray(
"RankArray"));
562 if(vertexRankArray !=
nullptr) {
563 for(
int i = 0; i < input->GetNumberOfPoints(); i++) {
564 vertexRankArray[i] = oldToNewRanks[vertexRankArray[i]];
568 = ttkUtils::GetPointer<int>(input->GetCellData()->GetArray(
"RankArray"));
569 if(cellRankArray !=
nullptr) {
570 for(
int i = 0; i < input->GetNumberOfCells(); i++) {
571 cellRankArray[i] = oldToNewRanks[cellRankArray[i]];
581 unsigned char *ghost,
584 if(rankArray !=
nullptr) {
585#ifdef TTK_ENABLE_OPENMP
586#pragma omp parallel for reduction(+ : ghostNumber)
594 if(ghost !=
nullptr) {
595#ifdef TTK_ENABLE_OPENMP
596#pragma omp parallel for reduction(+ : ghostNumber)
607 auto minmax = std::minmax_element(globalIds, globalIds + simplexNumber);
613 MPI_Allreduce(&realSimplexNumber, &globalSimplexNumber, 1,
614 ttk::getMPIType(realSimplexNumber), MPI_SUM, ttk::MPIcomm_);
616 &min, &globalMin, 1, ttk::getMPIType(min), MPI_MIN, ttk::MPIcomm_);
618 &max, &globalMax, 1, ttk::getMPIType(max), MPI_MAX, ttk::MPIcomm_);
620 return (globalSimplexNumber == globalMax + 1 && globalMin == 0);
625 std::unordered_map<ttk::SimplexId, ttk::SimplexId> &vertGtoL,
626 std::vector<int> &neighborRanks,
627 std::map<int, int> &neighborsToId) {
631 vtkNew<vtkIdTypeArray> vtkVertexIdentifiers{};
633 vtkNew<vtkIdTypeArray> vtkCellIdentifiers{};
634 vtkVertexIdentifiers->SetName(
"GlobalPointIds");
635 vtkVertexIdentifiers->SetNumberOfComponents(1);
636 vtkVertexIdentifiers->SetNumberOfTuples(input->GetNumberOfPoints());
637 vtkVertexIdentifiers->Fill(-1);
640 ttkUtils::GetPointer<ttk::LongSimplexId>(vtkVertexIdentifiers));
642 vtkCellIdentifiers->SetName(
"GlobalCellIds");
643 vtkCellIdentifiers->SetNumberOfComponents(1);
644 vtkCellIdentifiers->SetNumberOfTuples(input->GetNumberOfCells());
645 vtkCellIdentifiers->Fill(-1);
648 ttkUtils::GetPointer<ttk::LongSimplexId>(vtkCellIdentifiers));
650 int vertexNumber = input->GetNumberOfPoints();
652 int cellNumber = input->GetNumberOfCells();
656 double *boundingBox = input->GetBounds();
657 identifiers.setBounds(boundingBox);
658 identifiers.initializeNeighbors(boundingBox, neighborRanks, neighborsToId);
659 if(ttk::isRunningWithMPI()) {
660 switch(input->GetDataObjectType()) {
661 case VTK_UNSTRUCTURED_GRID:
662 case VTK_POLY_DATA: {
664 identifiers.setOutdatedGlobalPointIds(
665 ttkUtils::GetPointer<ttk::LongSimplexId>(
666 input->GetPointData()->GetGlobalIds()));
667 identifiers.setOutdatedGlobalCellIds(
668 ttkUtils::GetPointer<ttk::LongSimplexId>(
669 input->GetCellData()->GetGlobalIds()));
670 identifiers.setVertexRankArray(ttkUtils::GetPointer<int>(
671 input->GetPointData()->GetArray(
"RankArray")));
672 int *cellRankArray = ttkUtils::GetPointer<int>(
673 input->GetCellData()->GetArray(
"RankArray"));
674 identifiers.setCellRankArray(cellRankArray);
675 identifiers.setVertGhost(ttkUtils::GetPointer<unsigned char>(
676 input->GetPointData()->GetArray(
"vtkGhostType")));
677 unsigned char *cellGhost = ttkUtils::GetPointer<unsigned char>(
678 input->GetCellData()->GetArray(
"vtkGhostType"));
679 identifiers.setCellGhost(cellGhost);
680 vtkPointSet *pointSet = vtkPointSet::SafeDownCast(input);
681 identifiers.setPointSet(
static_cast<float *
>(
683 vtkCellArray *cells =
nullptr;
684 switch(input->GetDataObjectType()) {
685 case VTK_UNSTRUCTURED_GRID: {
686 auto dataSetAsUG = vtkUnstructuredGrid::SafeDownCast(input);
687 cells = dataSetAsUG->GetCells();
690 case VTK_POLY_DATA: {
691 auto dataSetAsPD = vtkPolyData::SafeDownCast(input);
693 = dataSetAsPD->GetNumberOfPolys() > 0 ? dataSetAsPD->GetPolys()
694 : dataSetAsPD->GetNumberOfLines() > 0 ? dataSetAsPD->GetLines()
695 : dataSetAsPD->GetVerts();
699 this->
printErr(
"Unable to get cells for `"
700 + std::string(input->GetClassName()) +
"`");
703 if(cells ==
nullptr) {
706 if(!cells->IsStorage64Bit()) {
707 if(cells->CanConvertTo64BitStorage()) {
708 this->
printWrn(
"Converting the cell array to 64-bit storage");
709 bool success = cells->ConvertTo64BitStorage();
712 "Error converting the provided cell array to 64-bit storage");
717 "Cannot convert the provided cell array to 64-bit storage");
722 identifiers.setConnectivity(ttkUtils::GetPointer<ttk::LongSimplexId>(
723 cells->GetConnectivityArray()));
725 std::vector<std::vector<ttk::SimplexId>> pointsToCells(vertexNumber);
726 vtkIdList *cellList = vtkIdList::New();
727 if(cellRankArray !=
nullptr) {
729 input->GetPointCells(i, cellList);
730 for(
int j = 0; j < cellList->GetNumberOfIds(); j++) {
732 pointsToCells[i].push_back(cellList->GetId(j));
738 input->GetPointCells(i, cellList);
739 for(
int j = 0; j < cellList->GetNumberOfIds(); j++) {
740 if(cellGhost[cellList->GetId(j)] == 0) {
741 pointsToCells[i].push_back(cellList->GetId(j));
746 identifiers.setPointsToCells(pointsToCells);
748 identifiers.initializeMPITypes();
749 identifiers.setVertGtoL(&vertGtoL);
750 vtkIdList *pointCell = vtkIdList::New();
751 input->GetCellPoints(0, pointCell);
752 int nbPoints = pointCell->GetNumberOfIds();
753 identifiers.setDomainDimension(nbPoints - 1);
754 identifiers.buildKDTree();
755 status = identifiers.executePolyData();
758 case VTK_IMAGE_DATA: {
759 vtkImageData *data = vtkImageData::SafeDownCast(input);
760 identifiers.setDims(data->GetDimensions());
761 identifiers.setSpacing(data->GetSpacing());
762 status = identifiers.executeImageData();
766 this->
printErr(
"Unable to triangulate `"
767 + std::string(input->GetClassName()) +
"`");
775 printErr(
"Global identifier generation failed");
781 input->GetPointData()->SetGlobalIds(vtkVertexIdentifiers);
782 input->GetCellData()->SetGlobalIds(vtkCellIdentifiers);
788 vtkNew<vtkGhostCellsGenerator> generator;
789 if(ttk::isRunningWithMPI()
790 && (!input->HasAnyGhostCells()
791 && ((input->GetPointData()->GetArray(
"RankArray") ==
nullptr)
792 || (input->GetCellData()->GetArray(
"RankArray") ==
nullptr)))) {
793 generator->SetInputData(input);
794 generator->BuildIfRequiredOff();
795 generator->SetNumberOfGhostLayers(1);
797 input->ShallowCopy(generator->GetOutputDataObject(0));
798 input->GetPointData()->AddArray(
799 generator->GetOutputDataObject(0)->GetGhostArray(0));
800 input->GetCellData()->AddArray(
801 generator->GetOutputDataObject(0)->GetGhostArray(1));
807 std::vector<int> &neighbors,
808 std::map<int, int> &neighToId,
814 if((input->GetDataObjectType() == VTK_POLY_DATA
815 || input->GetDataObjectType() == VTK_UNSTRUCTURED_GRID)) {
817 if((ttk::hasInitializedMPI()) && (ttk::isRunningWithMPI())) {
819 printWrn(
"The distribution by VTK of Unstructured");
820 printWrn(
"Grids and Poly Data has been reported");
821 printWrn(
"to be affected by bugs (at least up");
824 if((input->GetCellData()->GetGlobalIds() ==
nullptr)
825 || (input->GetPointData()->GetGlobalIds() ==
nullptr)) {
827 printWrn(
"=> Global identifiers may be incorrect.");
829 if((input->GetPointData()->GetArray(
"RankArray") ==
nullptr)
830 || (input->GetCellData()->GetArray(
"RankArray") ==
nullptr)) {
831 printWrn(
"=> Rank arrays may be incorrect.");
838 std::vector<int> &neighborRanks{
839 triangulation !=
nullptr ? triangulation->getNeighborRanks() : neighbors};
840 std::map<int, int> &neighborsToId{
841 triangulation !=
nullptr ? triangulation->getNeighborsToId() : neighToId};
843 double *boundingBox = input->GetBounds();
844 if(triangulation !=
nullptr) {
845 triangulation->createMetaGrid(boundingBox);
848 if(neighborRanks.empty()) {
849 ttk::preconditionNeighborsUsingBoundingBox(
850 boundingBox, neighborRanks, neighborsToId);
855 input->GetPointData()->GetGlobalIds());
857 input->GetCellData()->GetGlobalIds());
859 bool pointValidity{
false};
860 bool cellValidity{
false};
861 if((triangulation !=
nullptr
864 || triangulation ==
nullptr) {
865 if(globalPointIds !=
nullptr) {
866 unsigned char *ghostPoints = ttkUtils::GetPointer<unsigned char>(
867 input->GetPointData()->GetArray(
"vtkGhostType"));
868 int *vertexRankArray = ttkUtils::GetPointer<int>(
869 input->GetPointData()->GetArray(
"RankArray"));
871 globalPointIds, vertexNumber, ghostPoints, vertexRankArray);
873 if(pointValidity && globalCellIds !=
nullptr) {
875 unsigned char *ghostCells = ttkUtils::GetPointer<unsigned char>(
876 input->GetCellData()->GetArray(
"vtkGhostType"));
877 int *cellRankArray = ttkUtils::GetPointer<int>(
878 input->GetCellData()->GetArray(
"RankArray"));
880 globalCellIds, cellNumber, ghostCells, cellRankArray);
883 pointValidity =
true;
888 if(!pointValidity || !cellValidity) {
889 if(triangulation !=
nullptr) {
893 neighborRanks, neighborsToId);
896 std::unordered_map<ttk::SimplexId, ttk::SimplexId> vertGtoL{};
905 const auto pd{input->GetPointData()};
907 triangulation->
printWrn(
"No point data on input object");
911 triangulation->setVertsGlobalIds(
912 ttkUtils::GetPointer<ttk::LongSimplexId>(pd->GetGlobalIds()));
913 triangulation->setVertexGhostArray(
914 ttkUtils::GetPointer<unsigned char>(pd->GetArray(
"vtkGhostType")));
915 int *vertexRankArray = ttkUtils::GetPointer<int>(pd->GetArray(
"RankArray"));
916 if(vertexRankArray !=
nullptr) {
917 triangulation->setVertexRankArray(vertexRankArray);
919 triangulation->preconditionDistributedVertices();
922 const auto cd{input->GetCellData()};
924 triangulation->
printWrn(
"No cell data on input object");
928 triangulation->setCellsGlobalIds(
929 ttkUtils::GetPointer<ttk::LongSimplexId>(cd->GetGlobalIds()));
930 triangulation->setCellGhostArray(
931 ttkUtils::GetPointer<unsigned char>(cd->GetArray(
"vtkGhostType")));
932 int *cellRankArray = ttkUtils::GetPointer<int>(cd->GetArray(
"RankArray"));
933 if(cellRankArray !=
nullptr) {
934 triangulation->setCellRankArray(cellRankArray);
936 triangulation->preconditionDistributedCells();
944 vtkInformationVector **inputVector,
945 vtkInformationVector *outputVector) {
947 if(request->Has(vtkCompositeDataPipeline::REQUEST_DATA_OBJECT())) {
954 if(request->Has(vtkCompositeDataPipeline::REQUEST_INFORMATION())) {
961 if(request->Has(vtkCompositeDataPipeline::REQUEST_UPDATE_TIME())) {
969 vtkCompositeDataPipeline::REQUEST_TIME_DEPENDENT_INFORMATION())) {
970 this->
printMsg(
"Processing REQUEST_TIME_DEPENDENT_INFORMATION",
973 request, inputVector, outputVector);
977 if(request->Has(vtkCompositeDataPipeline::REQUEST_UPDATE_EXTENT())) {
984 if(request->Has(vtkCompositeDataPipeline::REQUEST_DATA_NOT_GENERATED())) {
991 if(request->Has(vtkCompositeDataPipeline::REQUEST_DATA())) {
995 if(ttk::hasInitializedMPI() && inputVector !=
nullptr) {
996 if(this->updateMPICommunicator(vtkDataSet::GetData(inputVector[0], 0))) {
1001 return this->
RequestData(request, inputVector, outputVector);
1004 this->
printErr(
"Unsupported pipeline pass:");
1005 request->Print(cout);
#define TTK_FORCE_USE(x)
Force the compiler to use the function/method parameter.
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
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))
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)
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)
int preconditionTriangulation(AbstractTriangulation *triangulation)
void setGlobalOrder(bool order)
int processScalarArray(const triangulationType *triangulation, ttk::SimplexId *orderArray, const DT *scalarArray, const size_t nVerts) const
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
std::string debugMsgNamePrefix_
void setDebugMsgPrefix(const std::string &prefix)
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
void setCellIdentifiers(ttk::LongSimplexId *cellIdentifiers)
void setVertexNumber(const SimplexId &vertexNumber)
void setCellNumber(const SimplexId &cellNumber)
int executeSequential()
Generates global ids for all data set type in sequential.
void setVertexIdentifiers(ttk::LongSimplexId *vertexIdentifiers)
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
Triangulation::Type getType() const
COMMON_EXPORTS int MPIsize_
COMMON_EXPORTS int MPIrank_
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.
int SimplexId
Identifier type for simplices of any dimension.
vtkStandardNewMacro(ttkAlgorithm)
vtkInformationKeyMacro(ttkAlgorithm, SAME_DATA_TYPE_AS_INPUT_PORT, Integer)
int prepOutput(vtkInformation *info, const std::string &className)
#define ttkTypeMacroAT(group0, group1, call)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)