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
|
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#include "vtkStructuredGrid.h"
#include "vtkCellData.h"
#include "vtkConstantArray.h"
#include "vtkGenericCell.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkMath.h"
#include "vtkObjectFactory.h"
#include "vtkPointData.h"
#include "vtkStructuredCellArray.h"
#include "vtkUnsignedCharArray.h"
VTK_ABI_NAMESPACE_BEGIN
vtkStandardNewMacro(vtkStructuredGrid);
vtkStandardExtendedNewMacro(vtkStructuredGrid);
//------------------------------------------------------------------------------
vtkStructuredGrid::vtkStructuredGrid()
{
this->DataDescription = VTK_EMPTY;
for (int idx = 0; idx < 3; ++idx)
{
this->Dimensions[idx] = 0;
}
const int extent[6] = { 0, -1, 0, -1, 0, -1 };
memcpy(this->Extent, extent, sizeof(extent));
this->Information->Set(vtkDataObject::DATA_EXTENT_TYPE(), VTK_3D_EXTENT);
this->Information->Set(vtkDataObject::DATA_EXTENT(), this->Extent, 6);
}
//------------------------------------------------------------------------------
vtkStructuredGrid::~vtkStructuredGrid() = default;
//------------------------------------------------------------------------------
// Copy the geometric and topological structure of an input structured grid.
void vtkStructuredGrid::CopyStructure(vtkDataSet* ds)
{
vtkStructuredGrid* sg = static_cast<vtkStructuredGrid*>(ds);
this->Superclass::CopyStructure(ds);
// set extent sets, extent, dimensions, and data description
this->SetExtent(sg->Extent);
if (ds->HasAnyBlankPoints())
{
// there is blanking
this->GetPointData()->AddArray(ds->GetPointGhostArray());
}
if (ds->HasAnyBlankCells())
{
// there is blanking
this->GetCellData()->AddArray(ds->GetCellGhostArray());
}
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::Initialize()
{
this->Superclass::Initialize();
if (this->Information)
{
this->SetDimensions(0, 0, 0);
}
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::BuildCells()
{
this->StructuredCells = vtkStructuredData::GetCellArray(this->Extent, false);
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::BuildCellTypes()
{
this->StructuredCellTypes = vtkStructuredData::GetCellTypesArray(this->Extent, false);
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::BuildImplicitStructures()
{
this->BuildCells();
this->BuildCellTypes();
}
//------------------------------------------------------------------------------
vtkCell* vtkStructuredGrid::GetCell(vtkIdType cellId)
{
this->GetCell(cellId, this->GenericCell);
return this->GenericCell->GetRepresentativeCell();
}
//------------------------------------------------------------------------------
vtkCell* vtkStructuredGrid::GetCell(int i, int j, int k)
{
int ijkMin[3] = { i, j, k };
const auto cellId = vtkStructuredData::ComputeCellId(this->Dimensions, ijkMin);
return this->GetCell(cellId);
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::GetCell(vtkIdType cellId, vtkGenericCell* cell)
{
// Make sure data is defined
if (!this->Points)
{
vtkErrorMacro(<< "No data");
return;
}
// see whether the cell is blanked
if (!this->IsCellVisible(cellId))
{
cell->SetCellTypeToEmptyCell();
return;
}
cell->SetCellType(this->StructuredCellTypes->GetValue(cellId));
this->StructuredCells->GetCellAtId(cellId, cell->PointIds);
this->Points->GetPoints(cell->PointIds, cell->Points);
}
//------------------------------------------------------------------------------
// Fast implementation of GetCellBounds(). Bounds are calculated without
// constructing a cell.
void vtkStructuredGrid::GetCellBounds(vtkIdType cellId, double bounds[6])
{
vtkIdType npts, pts[8];
this->StructuredCells->GetCellAtId(cellId, npts, pts);
vtkBoundingBox::ComputeBounds(this->Points, pts, npts, bounds);
}
//------------------------------------------------------------------------------
int vtkStructuredGrid::GetCellType(vtkIdType cellId)
{
// see whether the cell is blanked
return this->IsCellVisible(cellId) ? this->StructuredCellTypes->GetValue(cellId) : VTK_EMPTY_CELL;
}
//------------------------------------------------------------------------------
vtkIdType vtkStructuredGrid::GetCellSize(vtkIdType cellId)
{
// see whether the cell is blanked
return this->IsCellVisible(cellId) ? this->StructuredCells->GetCellSize(cellId) : 0;
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts,
vtkIdList* ptIds) VTK_SIZEHINT(pts, npts)
{
this->StructuredCells->GetCellAtId(cellId, npts, pts, ptIds);
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::GetCellPoints(vtkIdType cellId, vtkIdList* ptIds)
{
this->StructuredCells->GetCellAtId(cellId, ptIds);
}
namespace
{
class CellVisibility
{
public:
CellVisibility(vtkStructuredGrid* input)
: Input(input)
{
}
bool operator()(const vtkIdType id) { return !Input->IsCellVisible(id); }
private:
vtkStructuredGrid* Input;
};
} // anonymous namespace
//------------------------------------------------------------------------------
void vtkStructuredGrid::GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds)
{
int numPtIds = ptIds->GetNumberOfIds();
// Use special methods for speed
switch (numPtIds)
{
case 0:
cellIds->Reset();
return;
case 1:
case 2:
case 4: // vertex, edge, face neighbors
vtkStructuredData::GetCellNeighbors(cellId, ptIds, cellIds, this->Dimensions);
break;
default:
this->Superclass::GetCellNeighbors(cellId, ptIds, cellIds);
}
// If blanking, remove blanked cells.
if (this->GetPointGhostArray() || this->GetCellGhostArray())
{
vtkIdType* pCellIds = cellIds->GetPointer(0);
vtkIdType* end =
std::remove_if(pCellIds, pCellIds + cellIds->GetNumberOfIds(), CellVisibility(this));
cellIds->Resize(std::distance(pCellIds, end));
}
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::GetCellNeighbors(
vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int* seedLoc)
{
int numPtIds = ptIds->GetNumberOfIds();
// Use special methods for speed
switch (numPtIds)
{
case 0:
cellIds->Reset();
return;
case 1:
case 2:
case 4: // vertex, edge, face neighbors
vtkStructuredData::GetCellNeighbors(cellId, ptIds, cellIds, this->Dimensions, seedLoc);
break;
default:
this->Superclass::GetCellNeighbors(cellId, ptIds, cellIds);
}
// If blanking, remove blanked cells.
if (this->GetPointGhostArray() || this->GetCellGhostArray())
{
vtkIdType* pCellIds = cellIds->GetPointer(0);
vtkIdType* end =
std::remove_if(pCellIds, pCellIds + cellIds->GetNumberOfIds(), CellVisibility(this));
cellIds->Resize(std::distance(pCellIds, end));
}
}
//------------------------------------------------------------------------------
vtkStructuredCellArray* vtkStructuredGrid::GetCells()
{
return this->StructuredCells;
}
//------------------------------------------------------------------------------
vtkConstantArray<int>* vtkStructuredGrid::GetCellTypesArray()
{
return this->StructuredCellTypes;
}
//------------------------------------------------------------------------------
// Turn off a particular data point.
void vtkStructuredGrid::BlankPoint(vtkIdType ptId)
{
vtkUnsignedCharArray* ghosts = this->GetPointGhostArray();
if (!ghosts)
{
ghosts = this->AllocatePointGhostArray();
}
ghosts->SetValue(ptId, ghosts->GetValue(ptId) | vtkDataSetAttributes::HIDDENPOINT);
assert(!this->IsPointVisible(ptId));
}
//------------------------------------------------------------------------------
// Turn on a particular data point.
void vtkStructuredGrid::UnBlankPoint(vtkIdType ptId)
{
vtkUnsignedCharArray* ghosts = this->GetPointGhostArray();
if (ghosts)
{
ghosts->SetValue(ptId, ghosts->GetValue(ptId) & ~vtkDataSetAttributes::HIDDENPOINT);
}
assert(this->IsPointVisible(ptId));
}
//------------------------------------------------------------------------------
// Turn off a particular data cell.
void vtkStructuredGrid::BlankCell(vtkIdType cellId)
{
vtkUnsignedCharArray* ghosts = this->GetCellGhostArray();
if (!ghosts)
{
ghosts = this->AllocateCellGhostArray();
}
ghosts->SetValue(cellId, ghosts->GetValue(cellId) | vtkDataSetAttributes::HIDDENCELL);
assert(!this->IsCellVisible(cellId));
}
//------------------------------------------------------------------------------
// Turn on a particular data cell.
void vtkStructuredGrid::UnBlankCell(vtkIdType cellId)
{
vtkUnsignedCharArray* ghosts = this->GetCellGhostArray();
if (ghosts)
{
ghosts->SetValue(cellId, ghosts->GetValue(cellId) & ~vtkDataSetAttributes::HIDDENCELL);
}
}
//------------------------------------------------------------------------------
unsigned char vtkStructuredGrid::IsPointVisible(vtkIdType pointId)
{
return vtkStructuredData::IsPointVisible(pointId, this->GetPointGhostArray());
}
//------------------------------------------------------------------------------
// Return non-zero if the specified cell is visible (i.e., not blanked)
unsigned char vtkStructuredGrid::IsCellVisible(vtkIdType cellId)
{
return vtkStructuredData::IsCellVisible(cellId, this->Dimensions, this->DataDescription,
this->GetCellGhostArray(), this->GetPointGhostArray());
}
//------------------------------------------------------------------------------
bool vtkStructuredGrid::HasAnyBlankPoints()
{
return this->PointData->HasAnyGhostBitSet(vtkDataSetAttributes::HIDDENPOINT);
}
//------------------------------------------------------------------------------
bool vtkStructuredGrid::HasAnyBlankCells()
{
int cellBlanking = this->CellData->HasAnyGhostBitSet(vtkDataSetAttributes::HIDDENCELL);
return cellBlanking || this->HasAnyBlankPoints();
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::GetCellDims(int cellDims[3])
{
for (int i = 0; i < 3; ++i)
{
cellDims[i] = ((this->Dimensions[i] - 1) < 1) ? 1 : this->Dimensions[i] - 1;
}
}
//------------------------------------------------------------------------------
// Set dimensions of structured grid dataset.
void vtkStructuredGrid::SetDimensions(int i, int j, int k)
{
this->SetExtent(0, i - 1, 0, j - 1, 0, k - 1);
}
//------------------------------------------------------------------------------
// Set dimensions of structured grid dataset.
void vtkStructuredGrid::SetDimensions(const int dims[3])
{
this->SetExtent(0, dims[0] - 1, 0, dims[1] - 1, 0, dims[2] - 1);
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::SetExtent(VTK_FUTURE_CONST int extent[6])
{
int description = vtkStructuredData::SetExtent(extent, this->Extent);
if (description < 0) // improperly specified
{
vtkErrorMacro(<< "Bad Extent, retaining previous values");
}
if (description == VTK_UNCHANGED)
{
return;
}
vtkStructuredData::GetDimensionsFromExtent(extent, this->Dimensions);
this->DataDescription = description;
this->BuildImplicitStructures();
this->Modified();
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::SetExtent(int xMin, int xMax, int yMin, int yMax, int zMin, int zMax)
{
int extent[6];
extent[0] = xMin;
extent[1] = xMax;
extent[2] = yMin;
extent[3] = yMax;
extent[4] = zMin;
extent[5] = zMax;
this->SetExtent(extent);
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::GetDimensions(int dims[3])
{
const int* extent = this->Extent;
dims[0] = extent[1] - extent[0] + 1;
dims[1] = extent[3] - extent[2] + 1;
dims[2] = extent[5] - extent[4] + 1;
}
//------------------------------------------------------------------------------
unsigned long vtkStructuredGrid::GetActualMemorySize()
{
return this->Superclass::GetActualMemorySize();
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::ShallowCopy(vtkDataObject* dataObject)
{
vtkStructuredGrid* grid = vtkStructuredGrid::SafeDownCast(dataObject);
if (grid != nullptr)
{
this->InternalStructuredGridCopy(grid);
}
this->Superclass::ShallowCopy(dataObject);
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::DeepCopy(vtkDataObject* dataObject)
{
auto mkhold = vtkMemkindRAII(this->GetIsInMemkind());
vtkStructuredGrid* grid = vtkStructuredGrid::SafeDownCast(dataObject);
if (grid != nullptr)
{
this->InternalStructuredGridCopy(grid);
}
this->Superclass::DeepCopy(dataObject);
}
//------------------------------------------------------------------------------
// This copies all the local variables (but not objects).
void vtkStructuredGrid::InternalStructuredGridCopy(vtkStructuredGrid* src)
{
// set extent sets, extent, dimensions, and data description
this->SetExtent(src->Extent);
}
//------------------------------------------------------------------------------
// Override this method because of blanking
void vtkStructuredGrid::ComputeScalarRange()
{
if (this->GetMTime() > this->ScalarRangeComputeTime)
{
vtkDataArray* ptScalars = this->PointData->GetScalars();
vtkDataArray* cellScalars = this->CellData->GetScalars();
double ptRange[2];
double cellRange[2];
double s;
ptRange[0] = VTK_DOUBLE_MAX;
ptRange[1] = VTK_DOUBLE_MIN;
if (ptScalars)
{
vtkIdType num = this->GetNumberOfPoints();
for (vtkIdType id = 0; id < num; ++id)
{
if (this->IsPointVisible(id))
{
s = ptScalars->GetComponent(id, 0);
if (s < ptRange[0])
{
ptRange[0] = s;
}
if (s > ptRange[1])
{
ptRange[1] = s;
}
}
}
}
cellRange[0] = ptRange[0];
cellRange[1] = ptRange[1];
if (cellScalars)
{
vtkIdType num = this->GetNumberOfCells();
for (vtkIdType id = 0; id < num; ++id)
{
if (this->IsCellVisible(id))
{
s = cellScalars->GetComponent(id, 0);
if (s < cellRange[0])
{
cellRange[0] = s;
}
if (s > cellRange[1])
{
cellRange[1] = s;
}
}
}
}
this->ScalarRange[0] = (cellRange[0] >= VTK_DOUBLE_MAX ? 0.0 : cellRange[0]);
this->ScalarRange[1] = (cellRange[1] <= VTK_DOUBLE_MIN ? 1.0 : cellRange[1]);
this->ScalarRangeComputeTime.Modified();
}
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::Crop(const int* updateExtent)
{
// Do nothing for empty datasets:
for (int dim = 0; dim < 3; ++dim)
{
if (this->Extent[2 * dim] > this->Extent[2 * dim + 1])
{
vtkDebugMacro(<< "Refusing to crop empty dataset.");
return;
}
}
int i, j, k;
int uExt[6];
const int* extent = this->Extent;
// If the update extent is larger than the extent,
// we cannot do anything about it here.
for (i = 0; i < 3; ++i)
{
uExt[i * 2] = updateExtent[i * 2];
if (uExt[i * 2] < extent[i * 2])
{
uExt[i * 2] = extent[i * 2];
}
uExt[i * 2 + 1] = updateExtent[i * 2 + 1];
if (uExt[i * 2 + 1] > extent[i * 2 + 1])
{
uExt[i * 2 + 1] = extent[i * 2 + 1];
}
}
// If extents already match, then we need to do nothing.
if (extent[0] == uExt[0] && extent[1] == uExt[1] && extent[2] == uExt[2] &&
extent[3] == uExt[3] && extent[4] == uExt[4] && extent[5] == uExt[5])
{
return;
}
else
{
// Get the points. Protect against empty data objects.
vtkPoints* inPts = this->GetPoints();
if (inPts == nullptr)
{
return;
}
vtkDebugMacro(<< "Cropping Grid");
vtkStructuredGrid* newGrid = vtkStructuredGrid::New();
vtkPointData* inPD = this->GetPointData();
vtkCellData* inCD = this->GetCellData();
vtkPointData* outPD = newGrid->GetPointData();
vtkCellData* outCD = newGrid->GetCellData();
// Allocate necessary objects
//
newGrid->SetExtent(uExt);
int outSize = (uExt[1] - uExt[0] + 1) * (uExt[3] - uExt[2] + 1) * (uExt[5] - uExt[4] + 1);
vtkPoints* newPts = inPts->NewInstance();
newPts->SetDataType(inPts->GetDataType());
newPts->SetNumberOfPoints(outSize);
outPD->CopyAllocate(inPD, outSize, outSize);
outCD->CopyAllocate(inCD, outSize, outSize);
// Traverse this data and copy point attributes to output
vtkIdType newId = 0;
int inInc1 = (extent[1] - extent[0] + 1);
int inInc2 = inInc1 * (extent[3] - extent[2] + 1);
for (k = uExt[4]; k <= uExt[5]; ++k)
{
int kOffset = (k - extent[4]) * inInc2;
for (j = uExt[2]; j <= uExt[3]; ++j)
{
int jOffset = (j - extent[2]) * inInc1;
for (i = uExt[0]; i <= uExt[1]; ++i)
{
vtkIdType idx = (i - extent[0]) + jOffset + kOffset;
newPts->SetPoint(newId, inPts->GetPoint(idx));
outPD->CopyData(inPD, idx, newId++);
}
}
}
// Traverse input data and copy cell attributes to output
newId = 0;
inInc1 = (extent[1] - extent[0]);
inInc2 = inInc1 * (extent[3] - extent[2]);
for (k = uExt[4]; k < uExt[5]; ++k)
{
int kOffset = (k - extent[4]) * inInc2;
for (j = uExt[2]; j < uExt[3]; ++j)
{
int jOffset = (j - extent[2]) * inInc1;
for (i = uExt[0]; i < uExt[1]; ++i)
{
vtkIdType idx = (i - extent[0]) + jOffset + kOffset;
outCD->CopyData(inCD, idx, newId++);
}
}
}
this->SetExtent(uExt);
this->SetPoints(newPts);
newPts->Delete();
inPD->ShallowCopy(outPD);
inCD->ShallowCopy(outCD);
newGrid->Delete();
}
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
const int* dims = this->Dimensions;
os << indent << "Dimensions: (" << dims[0] << ", " << dims[1] << ", " << dims[2] << ")\n";
const int* extent = this->Extent;
os << indent << "Extent: " << extent[0] << ", " << extent[1] << ", " << extent[2] << ", "
<< extent[3] << ", " << extent[4] << ", " << extent[5] << endl;
os << ")\n";
}
//------------------------------------------------------------------------------
vtkStructuredGrid* vtkStructuredGrid::GetData(vtkInformation* info)
{
return info ? vtkStructuredGrid::SafeDownCast(info->Get(DATA_OBJECT())) : nullptr;
}
//------------------------------------------------------------------------------
vtkStructuredGrid* vtkStructuredGrid::GetData(vtkInformationVector* v, int i)
{
return vtkStructuredGrid::GetData(v->GetInformationObject(i));
}
//------------------------------------------------------------------------------
void vtkStructuredGrid::GetPoint(int i, int j, int k, double p[3], bool adjustForExtent)
{
if (i < this->Extent[0] || i > this->Extent[1] || j < this->Extent[2] || j > this->Extent[3] ||
k < this->Extent[4] || k > this->Extent[5])
{
vtkErrorMacro("ERROR: IJK coordinates are outside of grid extent!");
return; // out of bounds!
}
int pos[3];
pos[0] = i;
pos[1] = j;
pos[2] = k;
vtkIdType id;
if (adjustForExtent)
{
id = vtkStructuredData::ComputePointIdForExtent(this->Extent, pos);
}
else
{
id = vtkStructuredData::ComputePointId(this->Dimensions, pos);
}
this->GetPoint(id, p);
}
VTK_ABI_NAMESPACE_END
|