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 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892
|
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#include "vtkPConnectivityFilter.h"
#include "vtkAOSDataArrayTemplate.h"
#include "vtkArrayDispatch.h"
#include "vtkBoundingBox.h"
#include "vtkCellData.h"
#include "vtkCellIterator.h"
#include "vtkDataSet.h"
#include "vtkDataSetSurfaceFilter.h"
#include "vtkIdTypeArray.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkKdTreePointLocator.h"
#include "vtkMPICommunicator.h"
#include "vtkMPIController.h"
#include "vtkMath.h"
#include "vtkMultiProcessController.h"
#include "vtkNew.h"
#include "vtkObjectFactory.h"
#include "vtkPointData.h"
#include "vtkThreshold.h"
#include "vtkTypeList.h"
#include "vtkUnsignedCharArray.h"
#include "vtkUnstructuredGrid.h"
#include "vtkWeakPointer.h"
#include <map>
#include <numeric>
#include <set>
#include <vector>
VTK_ABI_NAMESPACE_BEGIN
namespace
{
typedef vtkTypeList::Unique<
vtkTypeList::Create<vtkAOSDataArrayTemplate<int>, vtkAOSDataArrayTemplate<unsigned long>,
vtkAOSDataArrayTemplate<char>, vtkAOSDataArrayTemplate<unsigned char>,
vtkAOSDataArrayTemplate<float>, vtkAOSDataArrayTemplate<double>>>::Result PointArrayTypes;
struct WorkerBase
{
WorkerBase(vtkMPIController* subController)
: SubController(subController)
{
}
/**
* MPI controller for ranks with data.
*/
vtkWeakPointer<vtkMPIController> SubController;
};
/**
* Worker for all ranks with data to exchange their bounding boxes.
*/
struct ExchangeBoundsWorker : public WorkerBase
{
ExchangeBoundsWorker(vtkMPIController* subController)
: WorkerBase(subController)
{
memset(this->Bounds, 0, sizeof(this->Bounds));
}
bool Execute(const double bounds[6], const vtkSmartPointer<vtkDataArray>& allBoundsArray)
{
memcpy(this->Bounds, bounds, 6 * sizeof(double));
using Dispatcher = vtkArrayDispatch::DispatchByArray<PointArrayTypes>;
return Dispatcher::Execute(allBoundsArray, *this);
}
template <class TArray>
void operator()(TArray* allBounds)
{
// Inflate the bounds a bit to deal with floating point precision.
double bounds[6];
vtkBoundingBox myBoundingBox(this->Bounds);
myBoundingBox.Inflate();
myBoundingBox.GetBounds(bounds);
typename TArray::ValueType typedBounds[6];
for (int i = 0; i < 6; ++i)
{
typedBounds[i] = static_cast<typename TArray::ValueType>(bounds[i]);
}
allBounds->SetNumberOfComponents(6);
allBounds->SetNumberOfTuples(this->SubController->GetNumberOfProcesses());
this->SubController->AllGather(
typedBounds, reinterpret_cast<typename TArray::ValueType*>(allBounds->GetVoidPointer(0)), 6);
}
protected:
// Input - Local data bounds
double Bounds[6];
};
/**
* Determine this rank's neighbors from the bounding box information.
*/
struct FindMyNeighborsWorker : public WorkerBase
{
FindMyNeighborsWorker(vtkMPIController* subController)
: WorkerBase(subController)
, MyNeighbors(nullptr)
{
memset(this->Bounds, 0, sizeof(this->Bounds));
}
bool Execute(const double bounds[6], const vtkSmartPointer<vtkDataArray>& allBoundsArray,
std::vector<int>& myNeighbors)
{
memcpy(this->Bounds, bounds, 6 * sizeof(double));
this->AllBoundsArray = allBoundsArray;
// Output
this->MyNeighbors = &myNeighbors;
this->DoExecute();
return true;
}
void DoExecute()
{
this->MyNeighbors->clear();
vtkBoundingBox myBoundingBox(this->Bounds);
myBoundingBox.Inflate();
double bounds[6];
myBoundingBox.GetBounds(bounds);
// Identify neighboring ranks.
int myRank = this->SubController->GetLocalProcessId();
for (int p = 0; p < this->SubController->GetNumberOfProcesses(); ++p)
{
if (p == myRank)
{
continue;
}
double potentialNeighborBounds[6];
this->AllBoundsArray->GetTuple(p, potentialNeighborBounds);
vtkBoundingBox potentialNeighborBoundingBox(potentialNeighborBounds);
if (myBoundingBox.Intersects(potentialNeighborBoundingBox))
{
this->MyNeighbors->push_back(p);
}
}
}
protected:
// Input - Local data bounds
double Bounds[6];
// Input - Bounds on all ranks
vtkWeakPointer<vtkDataArray> AllBoundsArray;
// Output - List of this rank's neighbors.
std::vector<int>* MyNeighbors;
};
/**
* Worker to gather up points and region ids to send to neighbors.
*/
struct AssemblePointsAndRegionIdsWorker : public WorkerBase
{
AssemblePointsAndRegionIdsWorker(vtkMPIController* subController)
: WorkerBase(subController)
, RegionStarts(nullptr)
, PointsForMyNeighbors(nullptr)
, RegionIdsForMyNeighbors(nullptr)
{
}
bool Execute(const std::vector<int>& regionStarts,
const vtkSmartPointer<vtkDataArray>& allBoundsArray,
const vtkSmartPointer<vtkPointSet>& localResult,
std::map<int, vtkSmartPointer<vtkDataArray>>& pointsForMyNeighbors,
std::map<int, vtkSmartPointer<vtkIdTypeArray>>& regionIdsForMyNeighbors)
{
this->RegionStarts = ®ionStarts;
this->LocalResult = localResult;
// Output
this->PointsForMyNeighbors = &pointsForMyNeighbors;
this->RegionIdsForMyNeighbors = ®ionIdsForMyNeighbors;
using Dispatcher = vtkArrayDispatch::DispatchByArray<PointArrayTypes>;
return Dispatcher::Execute(allBoundsArray, *this);
}
template <class TArray>
void operator()(TArray* allBounds)
{
// For all neighbors, gather up points and region IDs that they will
// potentially need. These are local points that fall within the bounding
// box of the neighbors.
this->PointsForMyNeighbors->clear();
this->RegionIdsForMyNeighbors->clear();
int myRank = this->SubController->GetLocalProcessId();
vtkPoints* outputPoints = this->LocalResult->GetPoints();
TArray* pointArray = TArray::SafeDownCast(outputPoints->GetData());
vtkPointData* outputPD = this->LocalResult->GetPointData();
vtkIdTypeArray* pointRegionIds = vtkIdTypeArray::SafeDownCast(outputPD->GetArray("RegionId"));
for (int p = 0; p < this->SubController->GetNumberOfProcesses(); ++p)
{
if (myRank == p)
{
continue;
}
TArray* typedPointsForMyNeighbor = TArray::New();
typedPointsForMyNeighbor->SetNumberOfComponents(3);
(*this->PointsForMyNeighbors)[p].TakeReference(typedPointsForMyNeighbor);
(*this->RegionIdsForMyNeighbors)[p] = vtkSmartPointer<vtkIdTypeArray>::New();
typename TArray::ValueType bb[6];
allBounds->GetTypedTuple(p, bb);
vtkBoundingBox neighborBB(bb[0], bb[1], bb[2], bb[3], bb[4], bb[5]);
for (vtkIdType id = 0; id < this->LocalResult->GetNumberOfPoints(); ++id)
{
typename TArray::ValueType pt[3];
pointArray->GetTypedTuple(id, pt);
double doublePt[3] = { static_cast<double>(pt[0]), static_cast<double>(pt[1]),
static_cast<double>(pt[2]) };
if (neighborBB.ContainsPoint(doublePt))
{
typedPointsForMyNeighbor->InsertNextTypedTuple(pt);
vtkIdType regionId =
pointRegionIds->GetTypedComponent(id, 0) + (*this->RegionStarts)[myRank];
(*this->RegionIdsForMyNeighbors)[p]->InsertNextTypedTuple(®ionId);
}
}
}
}
protected:
// Input - Starting index of the first region on each rank.
const std::vector<int>* RegionStarts;
// Input - Output from the local connectivity operation.
vtkWeakPointer<vtkPointSet> LocalResult;
// Output
std::map<int, vtkSmartPointer<vtkDataArray>>* PointsForMyNeighbors;
// Output
std::map<int, vtkSmartPointer<vtkIdTypeArray>>* RegionIdsForMyNeighbors;
};
/**
* Send and receive points to/from neighbors.
*/
struct SendReceivePointsWorker : public WorkerBase
{
SendReceivePointsWorker(vtkMPIController* subController)
: WorkerBase(subController)
, PointsFromMyNeighbors(nullptr)
, RegionIdsFromMyNeighbors(nullptr)
{
}
bool Execute(const vtkSmartPointer<vtkDataArray>& allBoundsArray,
const std::map<int, int>& sendLengths, const std::map<int, int>& recvLengths,
const std::map<int, vtkSmartPointer<vtkDataArray>>& pointsForMyNeighbors,
const std::map<int, vtkSmartPointer<vtkIdTypeArray>>& regionIdsForMyNeighbors,
std::map<int, vtkSmartPointer<vtkDataArray>>& pointsFromMyNeighbors,
std::map<int, vtkSmartPointer<vtkIdTypeArray>>& regionIdsFromMyNeighbors)
{
this->SendLengths = sendLengths;
this->RecvLengths = recvLengths;
this->PointsForMyNeighbors = pointsForMyNeighbors;
this->RegionIdsForMyNeighbors = regionIdsForMyNeighbors;
// Output
this->PointsFromMyNeighbors = &pointsFromMyNeighbors;
this->RegionIdsFromMyNeighbors = ®ionIdsFromMyNeighbors;
using Dispatcher = vtkArrayDispatch::DispatchByArray<PointArrayTypes>;
return Dispatcher::Execute(allBoundsArray, *this);
}
template <class TArray>
void operator()(TArray* vtkNotUsed(array))
{
const int PCF_POINTS_TAG = 194728;
const int PCF_REGIONIDS_TAG = 194729;
this->PointsFromMyNeighbors->clear();
this->RegionIdsFromMyNeighbors->clear();
std::map<int, vtkMPICommunicator::Request> sendRequestsPoints;
std::map<int, vtkMPICommunicator::Request> sendRequestsRegionIds;
std::vector<vtkMPICommunicator::Request> recvRequestsPoints(this->RecvLengths.size());
std::vector<vtkMPICommunicator::Request> recvRequestsRegionIds(this->RecvLengths.size());
// Receive neighbors' points.
int requestIdx = 0;
for (auto iter = this->RecvLengths.begin(); iter != this->RecvLengths.end(); ++iter)
{
int fromRank = iter->first;
int numFromRank = iter->second;
if (numFromRank > 0)
{
TArray* pfmn = TArray::New();
pfmn->SetNumberOfComponents(3);
pfmn->SetNumberOfTuples(numFromRank);
(*this->PointsFromMyNeighbors)[fromRank].TakeReference(pfmn);
this->SubController->NoBlockReceive(pfmn->GetPointer(0), 3 * numFromRank, fromRank,
PCF_POINTS_TAG, recvRequestsPoints[requestIdx]);
vtkIdTypeArray* idArray = vtkIdTypeArray::New();
idArray->SetNumberOfComponents(1);
idArray->SetNumberOfTuples(numFromRank);
(*this->RegionIdsFromMyNeighbors)[fromRank].TakeReference(idArray);
this->SubController->NoBlockReceive(idArray->GetPointer(0), numFromRank, fromRank,
PCF_REGIONIDS_TAG, recvRequestsRegionIds[requestIdx++]);
}
}
// Send points to neighbors
for (auto nbrIter = this->SendLengths.begin(); nbrIter != this->SendLengths.end(); ++nbrIter)
{
vtkIdType toRank = nbrIter->first;
vtkIdType numToRank = nbrIter->second;
if (numToRank > 0)
{
TArray* pfmn = TArray::SafeDownCast(this->PointsForMyNeighbors.at(toRank));
this->SubController->NoBlockSend(
pfmn->GetPointer(0), 3 * numToRank, toRank, PCF_POINTS_TAG, sendRequestsPoints[toRank]);
vtkIdTypeArray* idArray = this->RegionIdsForMyNeighbors.at(toRank);
this->SubController->NoBlockSend(idArray->GetPointer(0), numToRank, toRank,
PCF_REGIONIDS_TAG, sendRequestsRegionIds[toRank]);
}
}
this->SubController->WaitAll(requestIdx, recvRequestsPoints.data());
this->SubController->WaitAll(requestIdx, recvRequestsRegionIds.data());
}
protected:
// Input
std::map<int, int> SendLengths;
std::map<int, int> RecvLengths;
std::map<int, vtkSmartPointer<vtkDataArray>> PointsForMyNeighbors;
std::map<int, vtkSmartPointer<vtkIdTypeArray>> RegionIdsForMyNeighbors;
// Output
std::map<int, vtkSmartPointer<vtkDataArray>>* PointsFromMyNeighbors;
std::map<int, vtkSmartPointer<vtkIdTypeArray>>* RegionIdsFromMyNeighbors;
};
/**
* Exchange number of points going to each neighbor. No dispatch is needed for this function.
*/
void ExchangeNumberOfPointsToSend(vtkMPIController* subController,
const std::vector<int>& myNeighbors,
const std::map<int, vtkSmartPointer<vtkIdTypeArray>>& regionIdsForMyNeighbors,
std::map<int, int>& sendLengths, std::map<int, int>& recvLengths)
{
const int PCF_SIZE_EXCHANGE_TAG = 194727;
recvLengths.clear();
std::vector<vtkMPICommunicator::Request> recvRequests(myNeighbors.size());
int requestIdx = 0;
for (auto nbrIter = myNeighbors.begin(); nbrIter != myNeighbors.end(); ++nbrIter)
{
int fromRank = *nbrIter;
subController->NoBlockReceive(
&recvLengths[fromRank], 1, fromRank, PCF_SIZE_EXCHANGE_TAG, recvRequests[requestIdx++]);
}
std::map<int, vtkMPICommunicator::Request> sendRequests;
// Send number of points neighbors should expect to receive.
for (auto nbrIter = myNeighbors.begin(); nbrIter != myNeighbors.end(); ++nbrIter)
{
int toRank = *nbrIter;
sendLengths[toRank] = static_cast<int>(regionIdsForMyNeighbors.at(toRank)->GetNumberOfValues());
subController->NoBlockSend(
&sendLengths[toRank], 1, toRank, PCF_SIZE_EXCHANGE_TAG, sendRequests[toRank]);
}
subController->WaitAll(requestIdx, recvRequests.data());
}
} // end anonymous namespace
vtkStandardNewMacro(vtkPConnectivityFilter);
vtkPConnectivityFilter::vtkPConnectivityFilter() = default;
vtkPConnectivityFilter::~vtkPConnectivityFilter() = default;
int vtkPConnectivityFilter::RequestData(
vtkInformation* request, vtkInformationVector** inputVector, vtkInformationVector* outputVector)
{
// Get the input
vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
if (!inInfo)
{
return 0;
}
vtkDataSet* input = vtkDataSet::SafeDownCast(inInfo->Get(vtkDataObject::DATA_OBJECT()));
vtkMultiProcessController* globalController = vtkMultiProcessController::GetGlobalController();
// Check how many ranks have data. If it is only one, running the superclass
// RequestData is sufficient as no global data exchange and RegionId
// relabeling is needed. It is worth checking to avoid issuing an
// unnecessary warning when a dataset resides entirely on one process.
int numRanks = 1;
int myRank = 0;
int ranksWithCells = 0;
int hasCells = input->GetNumberOfCells() > 0 ? 1 : 0;
if (globalController)
{
globalController->AllReduce(&hasCells, &ranksWithCells, 1, vtkCommunicator::SUM_OP);
numRanks = globalController->GetNumberOfProcesses();
myRank = globalController->GetLocalProcessId();
}
// Compute local connectivity. If we are running in parallel, we need the full
// connectivity first, and will handle the extraction mode later.
int success = 1;
if (numRanks > 1 && ranksWithCells > 1)
{
if (this->ExtractionMode == VTK_EXTRACT_POINT_SEEDED_REGIONS ||
this->ExtractionMode == VTK_EXTRACT_CELL_SEEDED_REGIONS ||
this->ExtractionMode == VTK_EXTRACT_SPECIFIED_REGIONS)
{
vtkErrorMacro("ExtractionMode " << this->GetExtractionModeAsString()
<< " is not supported in " << this->GetClassName()
<< " when the number of ranks with data is greater than 1.");
return 1;
}
int saveScalarConnectivity = this->ScalarConnectivity;
int saveExtractionMode = this->ExtractionMode;
int saveColorRegions = this->ColorRegions;
int saveRegionIdAssignmentMode = this->RegionIdAssignmentMode;
bool compressArrays = this->GetCompressArrays();
// Overwrite custom member variables temporarily.
this->ScalarConnectivity = 0;
this->ExtractionMode = VTK_EXTRACT_ALL_REGIONS;
this->ColorRegions = 1;
this->RegionIdAssignmentMode = UNSPECIFIED;
this->CompressArraysOff();
// Invoke the connectivity algorithm in the superclass.
success = this->Superclass::RequestData(request, inputVector, outputVector);
this->ScalarConnectivity = saveScalarConnectivity;
this->ExtractionMode = saveExtractionMode;
this->ColorRegions = saveColorRegions;
this->RegionIdAssignmentMode = saveRegionIdAssignmentMode;
this->SetCompressArrays(compressArrays);
}
else
{
// Only 1 process, just invoke the superclass and return.
return this->Superclass::RequestData(request, inputVector, outputVector);
}
// Create a SubController.
vtkSmartPointer<vtkMPIController> subController;
subController.TakeReference(
vtkMPIController::SafeDownCast(globalController)->PartitionController(hasCells, 0));
// From here on we deal only with the SubController
numRanks = subController->GetNumberOfProcesses();
myRank = subController->GetLocalProcessId();
// Get the info objects
vtkInformation* outInfo = outputVector->GetInformationObject(0);
// Get the output
vtkPointSet* output = vtkPointSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
if (!output->GetPoints())
{
vtkErrorMacro("No points in data set");
success = 0;
}
// Check that all ranks succeeded in local connectivity.
int globalSuccess = 0;
subController->AllReduce(&success, &globalSuccess, 1, vtkCommunicator::MIN_OP);
if (!globalSuccess)
{
vtkErrorMacro("An error occurred on at least one process.");
return 0;
}
// Exchange number of regions. We assume the RegionIDs are contiguous.
int numRegions = this->GetNumberOfExtractedRegions();
std::vector<int> regionCounts(numRanks, 0);
std::vector<int> regionStarts(numRanks + 1, 0);
subController->AllGather(&numRegions, regionCounts.data(), 1);
// Compute starting region Ids on each rank
std::partial_sum(regionCounts.begin(), regionCounts.end(), regionStarts.begin() + 1);
vtkPointData* outputPD = output->GetPointData();
vtkIdTypeArray* pointRegionIds = vtkIdTypeArray::SafeDownCast(outputPD->GetArray("RegionId"));
// Exchange bounding boxes of the data on each rank. These are used to
// determine neighboring ranks and to minimize the number of points sent
// to neighboring processors.
vtkSmartPointer<vtkDataArray> allBoundsArray;
allBoundsArray.TakeReference(output->GetPoints()->GetData()->NewInstance());
ExchangeBoundsWorker exchangeBounds(subController);
if (!exchangeBounds.Execute(output->GetBounds(), allBoundsArray))
{
vtkErrorMacro("Unsupported points array type encountered when exchanging bounds.");
return 0;
}
// Identify neighboring ranks.
FindMyNeighborsWorker findMyNeighbors(subController);
std::vector<int> myNeighbors;
if (!findMyNeighbors.Execute(output->GetBounds(), allBoundsArray, myNeighbors))
{
vtkErrorMacro("Unsupported points array type encountered when finding neighbors.");
return 0;
}
AssemblePointsAndRegionIdsWorker assemblePointsAndRegionIds(subController);
std::map<int, vtkSmartPointer<vtkDataArray>> pointsForMyNeighbors;
std::map<int, vtkSmartPointer<vtkIdTypeArray>> regionIdsForMyNeighbors;
if (!assemblePointsAndRegionIds.Execute(
regionStarts, allBoundsArray, output, pointsForMyNeighbors, regionIdsForMyNeighbors))
{
vtkErrorMacro(
"Unsupported points array type encountered when assembling points and region ids.");
return 0;
}
std::map<int, int> sendLengths;
std::map<int, int> recvLengths;
ExchangeNumberOfPointsToSend(
subController, myNeighbors, regionIdsForMyNeighbors, sendLengths, recvLengths);
SendReceivePointsWorker sendReceivePoints(subController);
std::map<int, vtkSmartPointer<vtkDataArray>> pointsFromMyNeighbors;
std::map<int, vtkSmartPointer<vtkIdTypeArray>> regionIdsFromMyNeighbors;
if (!sendReceivePoints.Execute(allBoundsArray, sendLengths, recvLengths, pointsForMyNeighbors,
regionIdsForMyNeighbors, pointsFromMyNeighbors, regionIdsFromMyNeighbors))
{
vtkErrorMacro("Unsupported points array type encountered when sending and receiving points.");
return 0;
}
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
// Links from local region ids to remote region ids. Vector index is local
// region id, and the set contains linked remote ids.
typedef std::vector<std::set<vtkIdType>> RegionLinksType;
RegionLinksType links(regionStarts[numRanks]);
if (output->GetNumberOfPoints() > 0)
{
// Now resolve the points from our neighbors to local points if possible.
vtkNew<vtkKdTreePointLocator> localPointLocator;
localPointLocator->SetDataSet(output);
localPointLocator->BuildLocator();
// Map the local and remote ids
for (int rank = 0; rank < numRanks; ++rank)
{
if (rank == myRank || pointsFromMyNeighbors.count(rank) == 0)
{
continue;
}
for (vtkIdType ptId = 0; ptId < pointsFromMyNeighbors[rank]->GetNumberOfTuples(); ++ptId)
{
double x[3];
pointsFromMyNeighbors[rank]->GetTuple(ptId, x);
vtkIdType localId = localPointLocator->FindClosestPoint(x);
// Skip local ghost points as we do not need ghost-ghost links.
vtkUnsignedCharArray* pointGhostArray = output->GetPointGhostArray();
if (pointGhostArray &&
pointGhostArray->GetTypedComponent(localId, 0) & vtkDataSetAttributes::DUPLICATEPOINT)
{
continue;
}
double y[3];
output->GetPoints()->GetPoint(localId, y);
double dist2 = vtkMath::Distance2BetweenPoints(x, y);
if (dist2 > 1e-6)
{
// Nearest point is too far away, move on.
continue;
}
// Save association between local and remote ids
vtkIdTypeArray* localRegionIds =
vtkIdTypeArray::SafeDownCast(outputPD->GetArray("RegionId"));
vtkIdType localRegionId =
localRegionIds->GetTypedComponent(localId, 0) + regionStarts[myRank];
vtkIdType remoteRegionId = regionIdsFromMyNeighbors[rank]->GetTypedComponent(ptId, 0);
if (links[localRegionId].count(remoteRegionId) == 0)
{
// Link not already established. Add it here.
links[localRegionId].insert(remoteRegionId);
}
}
}
}
// Set up storage for gathering all links from all processors. This is an
// interleaved vector containing one regionId and its connected regionId.
std::vector<vtkIdType> localLinks;
for (size_t i = 0; i < links.size(); ++i)
{
for (std::set<vtkIdType>::iterator iter = std::begin(links[i]); iter != std::end(links[i]);
++iter)
{
localLinks.push_back(static_cast<vtkIdType>(i));
localLinks.push_back(*iter);
}
}
// Gather all the links on each rank. This is possibly suboptimal, but it
// avoids needing a connected components algorithm on a distributed graph.
vtkIdType localNumLinks = static_cast<vtkIdType>(localLinks.size());
std::vector<vtkIdType> linkCounts(numRanks, -1);
std::vector<vtkIdType> linkStarts(numRanks + 1, 0);
subController->AllGather(&localNumLinks, linkCounts.data(), 1);
// Compute starting region IDs on each rank
for (int i = 0; i < numRanks; ++i)
{
linkStarts[i + 1] = linkCounts[i] + linkStarts[i];
}
std::vector<vtkIdType> allLinks(linkStarts[numRanks]);
subController->AllGatherV(localLinks.data(), allLinks.data(),
static_cast<vtkIdType>(localLinks.size()), linkCounts.data(), linkStarts.data());
// Set up a graph of all the region-to-region links.
typedef struct _RegionNode
{
// Stored for relabeling step
vtkIdType OriginalRegionId;
// Current local region id
vtkIdType CurrentRegionId;
std::vector<vtkIdType> Links;
} RegionNode;
size_t linkIdx = 0;
std::vector<RegionNode> regionNodes(regionStarts[numRanks]);
for (vtkIdType regionId = 0; regionId < regionStarts[numRanks]; ++regionId)
{
regionNodes[regionId].OriginalRegionId = regionId;
regionNodes[regionId].CurrentRegionId = regionId;
while (linkIdx < allLinks.size() && allLinks[linkIdx] == regionId)
{
regionNodes[regionId].Links.push_back(allLinks[linkIdx + 1]);
linkIdx += 2;
}
}
// Now run connected components on this graph. The algorithm labels all
// connected nodes in the graph with the lowest region id in the connected
// component. This is a breadth-first algorithm. I'm not 100% sure that the
// multiple passes in the do-while loop are required, but I suspect there may
// be graph configurations where a single pass is not sufficient for the
// relabeling to converge.
bool componentChanged = false;
do
{
componentChanged = false;
for (std::vector<RegionNode>::iterator nodeIter = regionNodes.begin();
nodeIter != regionNodes.end(); ++nodeIter)
{
for (std::vector<vtkIdType>::iterator linkIter = nodeIter->Links.begin();
linkIter != nodeIter->Links.end(); ++linkIter)
{
vtkIdType linkedRegionId = *linkIter;
if (nodeIter->CurrentRegionId < regionNodes[linkedRegionId].CurrentRegionId)
{
regionNodes[linkedRegionId].CurrentRegionId = nodeIter->CurrentRegionId;
componentChanged = true;
}
}
}
} while (componentChanged);
// Collect all the current ids remaining after the connected components
// algorithm.
std::set<vtkIdType> currentRegionIds;
for (std::vector<RegionNode>::iterator nodeIter = regionNodes.begin();
nodeIter != regionNodes.end(); ++nodeIter)
{
currentRegionIds.insert(nodeIter->CurrentRegionId);
}
// Create a map from current region id after relabeling to a new, contiguous
// label. Maps current region id -> relabeled array.
std::map<vtkIdType, vtkIdType> relabeledRegionMap;
vtkIdType contiguousLabel = 0;
for (std::set<vtkIdType>::iterator setIter = currentRegionIds.begin();
setIter != currentRegionIds.end(); ++setIter)
{
relabeledRegionMap[*setIter] = contiguousLabel++;
}
// Now do the relabing to the contiguous region id.
std::vector<vtkIdType> regionIdMap(regionNodes.size(), -1);
for (std::vector<RegionNode>::iterator nodeIter = regionNodes.begin();
nodeIter != regionNodes.end(); ++nodeIter)
{
nodeIter->CurrentRegionId = relabeledRegionMap[nodeIter->CurrentRegionId];
}
// Relabel the points and cells according to the contiguous renumbering.
vtkCellData* outputCD = output->GetCellData();
vtkIdTypeArray* cellRegionIds = vtkIdTypeArray::SafeDownCast(outputCD->GetArray("RegionId"));
for (vtkIdType i = 0; i < output->GetNumberOfCells(); ++i)
{
// Offset the cellRegionId by the starting region id on this rank.
vtkIdType cellRegionId = cellRegionIds->GetValue(i) + regionStarts[myRank];
cellRegionIds->SetValue(i, regionNodes[cellRegionId].CurrentRegionId);
}
for (vtkIdType i = 0; i < output->GetNumberOfPoints(); ++i)
{
// Offset the pointRegionId by the starting region id on this rank.
vtkIdType pointRegionId = pointRegionIds->GetValue(i) + regionStarts[myRank];
pointRegionIds->SetValue(i, regionNodes[pointRegionId].CurrentRegionId);
}
// Sum up number of cells in each region.
vtkIdType numContiguousLabels = contiguousLabel;
std::vector<vtkIdType> localRegionSizes(numContiguousLabels, 0);
if (cellRegionIds)
{
// Iterate over cells and count how many are in different regions. Count only non-ghost cells.
vtkUnsignedCharArray* cellGhostArray = output->GetCellGhostArray();
for (vtkIdType i = 0; i < cellRegionIds->GetNumberOfValues(); ++i)
{
if (cellGhostArray &&
cellGhostArray->GetTypedComponent(i, 0) & vtkDataSetAttributes::DUPLICATECELL)
{
continue;
}
localRegionSizes[cellRegionIds->GetValue(i)]++;
}
}
// AllReduce to sum up the number of cells in each region on each process.
std::vector<vtkIdType> globalRegionSizes(numContiguousLabels, 0);
subController->AllReduce(localRegionSizes.data(), globalRegionSizes.data(), numContiguousLabels,
vtkCommunicator::SUM_OP);
// Store the region sizes
this->RegionSizes->Reset();
this->RegionSizes->SetNumberOfComponents(1);
this->RegionSizes->SetNumberOfTuples(numContiguousLabels);
for (vtkIdType i = 0; i < numContiguousLabels; ++i)
{
this->RegionSizes->SetTypedTuple(i, &globalRegionSizes[i]);
}
// Potentially reorder RegionIds in the output arrays.
this->OrderRegionIds(pointRegionIds, cellRegionIds);
if (this->ExtractionMode == VTK_EXTRACT_LARGEST_REGION ||
this->ExtractionMode == VTK_EXTRACT_CLOSEST_POINT_REGION)
{
double threshold = 0.0;
if (this->ExtractionMode == VTK_EXTRACT_LARGEST_REGION)
{
vtkIdType largestRegionCount = 0;
vtkIdType largestRegionId = 0;
for (vtkIdType i = 0; i < this->RegionSizes->GetNumberOfTuples(); ++i)
{
vtkIdType candidateCount = this->RegionSizes->GetValue(i);
if (candidateCount > largestRegionCount)
{
largestRegionCount = candidateCount;
largestRegionId = i;
}
}
threshold = largestRegionId;
}
else if (this->ExtractionMode == VTK_EXTRACT_CLOSEST_POINT_REGION)
{
// Find point closest to the desired point
double minDist2 = VTK_DOUBLE_MAX;
vtkIdType minId = 0;
for (vtkIdType i = 0; i < output->GetNumberOfPoints(); ++i)
{
double x[3];
output->GetPoint(i, x);
double dist2 = vtkMath::Distance2BetweenPoints(x, this->ClosestPoint);
if (dist2 < minDist2)
{
minDist2 = dist2;
minId = i;
}
}
// AllReduce to find the global minDist2.
double globalMinDist2 = VTK_DOUBLE_MAX;
subController->AllReduce(&minDist2, &globalMinDist2, 1, vtkCommunicator::MIN_OP);
int minDist2Rank = 0;
vtkIdType minDist2Region = 0;
if (fabs(minDist2 - globalMinDist2) < 1e-9)
{
minDist2Rank = myRank;
minDist2Region = pointRegionIds->GetValue(minId);
}
// Broadcast the rank of who has the minimum distance
int globalMinDist2Rank = 0;
subController->AllReduce(&minDist2Rank, &globalMinDist2Rank, 1, vtkCommunicator::MAX_OP);
// Get the id of the region nearest the point and use that in the
// threshold filter below.
subController->Broadcast(&minDist2Region, 1, globalMinDist2Rank);
threshold = minDist2Region;
}
// Now extract only the cells that have the desired id.
vtkNew<vtkThreshold> thresholder;
thresholder->SetInputData(output);
thresholder->SetThresholdFunction(vtkThreshold::THRESHOLD_BETWEEN);
thresholder->SetLowerThreshold(threshold);
thresholder->SetUpperThreshold(threshold);
thresholder->SetInputArrayToProcess(
0, 0, 0, vtkDataObject::FIELD_ASSOCIATION_CELLS, "RegionId");
thresholder->Update();
if (output->IsA("vtkPolyData"))
{
// It's too bad we have to do this, but vtkThreshold produces
// vtkUnstructuredGrid output.
vtkNew<vtkDataSetSurfaceFilter> surfaceFilter;
surfaceFilter->SetInputConnection(thresholder->GetOutputPort());
surfaceFilter->PassThroughCellIdsOff();
surfaceFilter->PassThroughPointIdsOff();
surfaceFilter->Update();
output->ShallowCopy(surfaceFilter->GetOutput());
}
else
{
// Output is an unstructured grid
output->DeepCopy(thresholder->GetOutput());
}
}
if (!this->ColorRegions)
{
// No coloring desired. Remove the RegionId arrays.
outputPD->RemoveArray("RegionId");
outputCD->RemoveArray("RegionId");
}
else
{
this->AddRegionsIds(output, pointRegionIds, cellRegionIds);
}
return 1;
}
void vtkPConnectivityFilter::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
}
VTK_ABI_NAMESPACE_END
|