VTK  9.4.20241217
vtkPolyData.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
165#ifndef vtkPolyData_h
166#define vtkPolyData_h
167
168#include "vtkCommonDataModelModule.h" // For export macro
169#include "vtkPointSet.h"
170#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
171
172#include "vtkCellArray.h" // Needed for inline methods
173#include "vtkCellLinks.h" // Needed for inline methods
174#include "vtkPolyDataInternals.h" // Needed for inline methods
175
176VTK_ABI_NAMESPACE_BEGIN
177struct vtkPolyDataDummyContainter;
179
180class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALAUTO vtkPolyData : public vtkPointSet
181{
182public:
183 static vtkPolyData* New();
185
186 vtkTypeMacro(vtkPolyData, vtkPointSet);
187 void PrintSelf(ostream& os, vtkIndent indent) override;
188
192 int GetDataObjectType() override { return VTK_POLY_DATA; }
193
197 void CopyStructure(vtkDataSet* ds) override;
198
200
203 vtkIdType GetNumberOfCells() override;
205 vtkCell* GetCell(vtkIdType cellId) override;
206 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
207 int GetCellType(vtkIdType cellId) override;
208 vtkIdType GetCellSize(vtkIdType cellId) override;
209 void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
210 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
212
220 void CopyCells(vtkPolyData* pd, vtkIdList* idList, vtkIncrementalPointLocator* locator = nullptr);
221
225 void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override;
226
231 void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override;
232
252
258 void GetCellsBounds(double bounds[6]);
259
266 void Squeeze() override;
267
271 int GetMaxCellSize() override;
272
274
281
288
293
299
304
310
315
321
326
333
335
338 vtkIdType GetNumberOfVerts() { return (this->Verts ? this->Verts->GetNumberOfCells() : 0); }
339 vtkIdType GetNumberOfLines() { return (this->Lines ? this->Lines->GetNumberOfCells() : 0); }
340 vtkIdType GetNumberOfPolys() { return (this->Polys ? this->Polys->GetNumberOfCells() : 0); }
341 vtkIdType GetNumberOfStrips() { return (this->Strips ? this->Strips->GetNumberOfCells() : 0); }
343
353 bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize);
354
364 bool AllocateEstimate(vtkIdType numVerts, vtkIdType maxVertSize, vtkIdType numLines,
365 vtkIdType maxLineSize, vtkIdType numPolys, vtkIdType maxPolySize, vtkIdType numStrips,
366 vtkIdType maxStripSize);
367
377 bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize);
378
389 bool AllocateExact(vtkIdType numVerts, vtkIdType vertConnSize, vtkIdType numLines,
390 vtkIdType lineConnSize, vtkIdType numPolys, vtkIdType polyConnSize, vtkIdType numStrips,
391 vtkIdType stripConnSize);
392
402
412 bool AllocateProportional(vtkPolyData* pd, double ratio);
413
420 void Allocate(vtkIdType numCells = 1000, int vtkNotUsed(extSize) = 1000)
421 {
422 this->AllocateExact(numCells, numCells);
423 }
424
435 void Allocate(vtkPolyData* inPolyData, vtkIdType numCells = 1000, int vtkNotUsed(extSize) = 1000)
436 {
437 this->AllocateProportional(
438 inPolyData, static_cast<double>(numCells) / inPolyData->GetNumberOfCells());
439 }
440
448 vtkIdType InsertNextCell(int type, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
449
458
463 void Reset();
464
473
477 bool NeedToBuildCells() { return this->Cells == nullptr; }
478
485 void BuildLinks(int initialSize = 0);
486
488
494
501
506
508
512 void GetPointCells(vtkIdType ptId, vtkIdType& ncells, vtkIdType*& cells)
513 VTK_SIZEHINT(cells, ncells);
515
522
534 unsigned char GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts)
535 VTK_SIZEHINT(pts, npts);
536
551 void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds)
552 VTK_SIZEHINT(pts, npts) override;
553
558 int IsTriangle(int v1, int v2, int v3);
559
568
573 int IsPointUsedByCell(vtkIdType ptId, vtkIdType cellId);
574
583 void ReplaceCell(vtkIdType cellId, vtkIdList* ids);
584 void ReplaceCell(vtkIdType cellId, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
588
597 void ReplaceCellPoint(vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId);
598 void ReplaceCellPoint(
599 vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId, vtkIdList* cellPointIds);
601
606 void ReverseCell(vtkIdType cellId);
607
609
613 void DeletePoint(vtkIdType ptId);
614 void DeleteCell(vtkIdType cellId);
616
626
628
637 vtkIdType InsertNextLinkedPoint(double x[3], int numLinks);
639
646 vtkIdType InsertNextLinkedCell(int type, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
647
657 void ReplaceLinkedCell(vtkIdType cellId, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
658
666 void RemoveCellReference(vtkIdType cellId);
667
675 void AddCellReference(vtkIdType cellId);
676
685
694
700 void ResizeCellList(vtkIdType ptId, int size);
701
705 void Initialize() override;
706
708
711 virtual int GetPiece();
712 virtual int GetNumberOfPieces();
714
718 virtual int GetGhostLevel();
719
728 unsigned long GetActualMemorySize() override;
729
731
734 void ShallowCopy(vtkDataObject* src) override;
735 void DeepCopy(vtkDataObject* src) override;
737
745
747
753
772 enum
773 {
774 ERR_NO_SUCH_FIELD = -4,
775 ERR_INCORRECT_FIELD = -3,
776 ERR_NON_MANIFOLD_STAR = -2,
777 REGULAR_POINT = -1,
778 MINIMUM = 0,
779 SADDLE = 1,
780 MAXIMUM = 2
781 };
782
784 int GetScalarFieldCriticalIndex(vtkIdType pointId, int fieldId);
785 int GetScalarFieldCriticalIndex(vtkIdType pointId, const char* fieldName);
786
795
800
810 unsigned char GetCell(vtkIdType cellId, const vtkIdType*& pts);
811
812protected:
814 ~vtkPolyData() override;
815
817
820
821 vtkCellArray* GetCellArrayInternal(TaggedCellId tag);
822
823 // points inherited
824 // point data (i.e., scalars, vectors, normals, tcoords) inherited
829
830 // supporting structures for more complex topological operations
831 // built only when necessary
834
836
837 // dummy static member below used as a trick to simplify traversal
838 static vtkPolyDataDummyContainter DummyContainer;
839
840 // Take into account only points that belong to at least one cell.
841 double CellsBounds[6];
842
844
845private:
846 void Cleanup();
847
848 vtkPolyData(const vtkPolyData&) = delete;
849 void operator=(const vtkPolyData&) = delete;
850};
851
852//------------------------------------------------------------------------------
854{
855 return (this->GetNumberOfVerts() + this->GetNumberOfLines() + this->GetNumberOfPolys() +
856 this->GetNumberOfStrips());
857}
858
859//------------------------------------------------------------------------------
861{
862 if (!this->Cells)
863 {
864 this->BuildCells();
865 }
866 return static_cast<int>(this->Cells->GetTag(cellId).GetCellType());
867}
868
869//------------------------------------------------------------------------------
871{
872 if (!this->Cells)
873 {
874 this->BuildCells();
875 }
876 switch (this->GetCellType(cellId))
877 {
878 case VTK_EMPTY_CELL:
879 return 0;
880 case VTK_VERTEX:
881 return 1;
882 case VTK_LINE:
883 return 2;
884 case VTK_TRIANGLE:
885 return 3;
886 case VTK_QUAD:
887 return 4;
888 case VTK_POLY_VERTEX:
889 return this->Verts ? this->Verts->GetCellSize(this->GetCellIdRelativeToCellArray(cellId)) : 0;
890 case VTK_POLY_LINE:
891 return this->Lines ? this->Lines->GetCellSize(this->GetCellIdRelativeToCellArray(cellId)) : 0;
892 case VTK_POLYGON:
893 return this->Polys ? this->Polys->GetCellSize(this->GetCellIdRelativeToCellArray(cellId)) : 0;
895 return this->Strips ? this->Strips->GetCellSize(this->GetCellIdRelativeToCellArray(cellId))
896 : 0;
897 }
898 vtkWarningMacro(<< "Cell type not supported.");
899 return 0;
900}
901
902//------------------------------------------------------------------------------
904{
905 vtkIdType npts;
906 const vtkIdType* pts;
907
908 this->GetCellPoints(cellId, npts, pts);
909 for (vtkIdType i = 0; i < npts; i++)
910 {
911 if (pts[i] == ptId)
912 {
913 return 1;
914 }
915 }
916
917 return 0;
918}
919
920//------------------------------------------------------------------------------
922{
923 static_cast<vtkCellLinks*>(this->Links.Get())->DeletePoint(ptId);
924}
925
926//------------------------------------------------------------------------------
928{
929 this->Cells->GetTag(cellId).MarkDeleted();
930}
931
932//------------------------------------------------------------------------------
934{
935 const vtkIdType* pts;
936 vtkIdType npts;
937
938 this->GetCellPoints(cellId, npts, pts);
939 auto links = static_cast<vtkCellLinks*>(this->Links.Get());
940 for (vtkIdType i = 0; i < npts; i++)
941 {
942 links->RemoveCellReference(cellId, pts[i]);
943 }
944}
945
946//------------------------------------------------------------------------------
948{
949 const vtkIdType* pts;
950 vtkIdType npts;
951
952 this->GetCellPoints(cellId, npts, pts);
953 auto links = static_cast<vtkCellLinks*>(this->Links.Get());
954 for (vtkIdType i = 0; i < npts; i++)
955 {
956 links->AddCellReference(cellId, pts[i]);
957 }
958}
959
960//------------------------------------------------------------------------------
961inline void vtkPolyData::ResizeCellList(vtkIdType ptId, int size)
962{
963 static_cast<vtkCellLinks*>(this->Links.Get())->ResizeCellList(ptId, size);
964}
965
966//------------------------------------------------------------------------------
968{
969 switch (tag.GetTarget())
970 {
972 return this->Verts;
974 return this->Lines;
976 return this->Polys;
978 return this->Strips;
979 }
980 return nullptr; // unreachable
981}
982
983//------------------------------------------------------------------------------
984inline void vtkPolyData::ReplaceCellPoint(vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId)
985{
987 this->ReplaceCellPoint(cellId, oldPtId, newPtId, ids);
988}
989
990//------------------------------------------------------------------------------
992 vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId, vtkIdList* cellPointIds)
993{
994 if (!this->Cells)
995 {
996 this->BuildCells();
997 }
998 vtkIdType npts;
999 const vtkIdType* pts;
1000 this->GetCellPoints(cellId, npts, pts, cellPointIds);
1001 for (vtkIdType i = 0; i < npts; i++)
1002 {
1003 if (pts[i] == oldPtId)
1004 {
1005 const TaggedCellId tag = this->Cells->GetTag(cellId);
1006 vtkCellArray* cells = this->GetCellArrayInternal(tag);
1007 cells->ReplaceCellPointAtId(tag.GetCellId(), i, newPtId);
1008 break;
1009 }
1010 }
1011}
1012
1013//------------------------------------------------------------------------------
1014inline unsigned char vtkPolyData::GetCellPoints(
1015 vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts)
1016{
1017 if (!this->Cells)
1018 {
1019 this->BuildCells();
1020 }
1021
1022 const TaggedCellId tag = this->Cells->GetTag(cellId);
1023 if (tag.IsDeleted())
1024 {
1025 npts = 0;
1026 pts = nullptr;
1027 return VTK_EMPTY_CELL;
1028 }
1029
1030 vtkCellArray* cells = this->GetCellArrayInternal(tag);
1031 cells->GetCellAtId(tag.GetCellId(), npts, pts);
1032 return tag.GetCellType();
1033}
1034
1035//------------------------------------------------------------------------------
1037 vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds)
1038{
1039 if (!this->Cells)
1040 {
1041 this->BuildCells();
1042 }
1043
1044 const TaggedCellId tag = this->Cells->GetTag(cellId);
1045 if (tag.IsDeleted())
1046 {
1047 npts = 0;
1048 pts = nullptr;
1049 }
1050 else
1051 {
1052 vtkCellArray* cells = this->GetCellArrayInternal(tag);
1053 cells->GetCellAtId(tag.GetCellId(), npts, pts, ptIds);
1054 }
1055}
1056
1057VTK_ABI_NAMESPACE_END
1058#endif
object to represent cell connectivity
vtkIdType GetCellSize(vtkIdType cellId) const override
Return the size of the cell at cellId.
void GetCellAtId(vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints, vtkIdList *ptIds) override
Return the point ids for the cell at cellId.
void ReplaceCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex, vtkIdType newPointId)
Replaces the pointId at cellPointIndex of a cell with newPointId.
abstract class to specify cell behavior
Definition vtkCell.h:130
abstract superclass for arrays of numeric data
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
virtual vtkCell * GetCell(vtkIdType cellId)=0
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
Detect and break reference loops.
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:133
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Allocate and hold a VTK object.
Definition vtkNew.h:167
concrete class for storing a set of points
Definition vtkPointSet.h:98
void GetCellPoints(vtkIdType, vtkIdList *idList) override
This method resets parameter idList, as there is no cell in a vtkPointSet.
vtkIdType GetCellSize(vtkIdType) override
This method always returns 1, as all cells are point in a pure vtkPointSet.
vtkIdType GetNumberOfCells() override
This method always returns 0, as there are no cells in a vtkPointSet.
int GetCellType(vtkIdType) override
This method always returns VTK_EMPTY_CELL, as there is no cell in a vtkPointSet.
concrete dataset represents vertices, lines, polygons, and triangle strips
vtkCellArray * GetStrips()
Get the cell array defining triangle strips.
vtkIdType InsertNextCell(int type, int npts, const vtkIdType pts[])
Insert a cell of type VTK_VERTEX, VTK_POLY_VERTEX, VTK_LINE, VTK_POLY_LINE, VTK_TRIANGLE,...
static vtkPolyData * ExtendedNew()
void GetCell(vtkIdType cellId, vtkGenericCell *cell) override
Standard vtkDataSet interface.
void Squeeze() override
Recover extra allocated memory when creating data whose initial size is unknown.
bool NeedToBuildCells()
Check if BuildCells is needed.
bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize)
Preallocate memory for the internal cell arrays.
int GetScalarFieldCriticalIndex(vtkIdType pointId, const char *fieldName)
vtkCell * GetCell(vtkIdType cellId) override
Standard vtkDataSet interface.
void SetPolys(vtkCellArray *p)
Set the cell array defining polygons.
vtkCellArray * GetCellArrayInternal(TaggedCellId tag)
vtkIdType GetCellIdRelativeToCellArray(vtkIdType cellId)
Maps the cell at position cellId inside the vtkPolyData to its location in the corresponding cell arr...
void GetCellsBounds(double bounds[6])
Get the cells bounds.
void RemoveCellReference(vtkIdType cellId)
Remove all references to cell in cell structure.
static vtkPolyData * GetData(vtkInformationVector *v, int i=0)
Retrieve an instance of this class from an information object.
void ComputeCellsBounds()
Compute the (X, Y, Z) bounds of the data.
vtkIdType GetNumberOfLines()
Return the number of primitives of a particular type held.
bool AllocateEstimate(vtkIdType numVerts, vtkIdType maxVertSize, vtkIdType numLines, vtkIdType maxLineSize, vtkIdType numPolys, vtkIdType maxPolySize, vtkIdType numStrips, vtkIdType maxStripSize)
Preallocate memory for the internal cell arrays.
void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds) override
Copy a cells point ids into list provided.
int IsTriangle(int v1, int v2, int v3)
Given three vertices, determine whether it's a triangle.
void SetLines(vtkCellArray *l)
Set the cell array defining lines.
void SetVerts(vtkCellArray *v)
Set the cell array defining vertices.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for type information and printing.
bool AllocateExact(vtkIdType numVerts, vtkIdType vertConnSize, vtkIdType numLines, vtkIdType lineConnSize, vtkIdType numPolys, vtkIdType polyConnSize, vtkIdType numStrips, vtkIdType stripConnSize)
Preallocate memory for the internal cell arrays.
void Initialize() override
Restore object to initial state.
void AddReferenceToCell(vtkIdType ptId, vtkIdType cellId)
Add a reference to a cell in a particular point's link list.
void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds) override
Standard vtkDataSet interface.
~vtkPolyData() override
int GetScalarFieldCriticalIndex(vtkIdType pointId, int fieldId)
void ReplaceLinkedCell(vtkIdType cellId, int npts, const vtkIdType pts[])
Replace one cell with another in cell structure.
vtkGetSmartPointerMacro(Links, vtkAbstractCellLinks)
Set/Get the links that were created possibly without using BuildLinks.
void CopyStructure(vtkDataSet *ds) override
Copy the geometric and topological structure of an input poly data object.
vtkNew< vtkIdList > LegacyBuffer
vtkSetSmartPointerMacro(Links, vtkAbstractCellLinks)
Set/Get the links that were created possibly without using BuildLinks.
int GetCellType(vtkIdType cellId) override
Standard vtkDataSet interface.
void GetCellBounds(vtkIdType cellId, double bounds[6]) override
Standard vtkDataSet interface.
void ReplaceCell(vtkIdType cellId, vtkIdList *ids)
Replace the points defining cell "cellId" with a new set of points.
virtual int GetGhostLevel()
Get the ghost level.
vtkMTimeType GetMTime() override
Get MTime which also considers its cell array MTime.
vtkSmartPointer< vtkCellArray > Verts
vtkIdType GetNumberOfStrips()
Return the number of primitives of a particular type held.
int GetScalarFieldCriticalIndex(vtkIdType pointId, vtkDataArray *scalarField)
bool AllocateCopy(vtkPolyData *pd)
Preallocate memory for the internal cell arrays such that they are the same size as those in pd.
vtkIdType InsertNextLinkedCell(int type, int npts, const vtkIdType pts[])
Add a new cell to the cell data structure (after cell pointers have been built).
vtkIdType InsertNextCell(int type, vtkIdList *pts)
Insert a cell of type VTK_VERTEX, VTK_POLY_VERTEX, VTK_LINE, VTK_POLY_LINE, VTK_TRIANGLE,...
vtkIdType GetNumberOfPolys()
Return the number of primitives of a particular type held.
void GetPointCells(vtkIdType ptId, vtkIdType &ncells, vtkIdType *&cells)
Special (efficient) operations on poly data.
vtkSmartPointer< vtkAbstractCellLinks > Links
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
vtkSmartPointer< vtkCellArray > Strips
void DeleteCell(vtkIdType cellId)
Mark a point/cell as deleted from this vtkPolyData.
vtkIdType InsertNextLinkedPoint(double x[3], int numLinks)
Add a point to the cell data structure (after cell pointers have been built).
int GetMaxCellSize() override
Return the maximum cell size in this poly data.
void RemoveGhostCells()
This method will remove any cell that is marked as ghost (has the vtkDataSetAttributes::DUPLICATECELL...
int GetMinSpatialDimension() override
Get the maximum/minimum spatial dimensionality of the data which is the maximum/minimum dimension of ...
void ResizeCellList(vtkIdType ptId, int size)
Resize the list of cells using a particular point.
virtual int GetPiece()
Get the piece and the number of pieces.
vtkTimeStamp CellsBoundsTime
int GetMaxSpatialDimension() override
Get the maximum/minimum spatial dimensionality of the data which is the maximum/minimum dimension of ...
vtkIdType GetNumberOfCells() override
Standard vtkDataSet interface.
bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize)
Preallocate memory for the internal cell arrays.
vtkIdType InsertNextLinkedPoint(int numLinks)
Add a point to the cell data structure (after cell pointers have been built).
unsigned char GetCell(vtkIdType cellId, const vtkIdType *&pts)
Get a pointer to the cell, ie [npts pid1 .
vtkCellArray * GetVerts()
Get the cell array defining vertices.
vtkSmartPointer< vtkCellArray > Polys
void Reset()
Begin inserting data all over again.
static vtkPolyData * New()
unsigned long GetActualMemorySize() override
Return the actual size of the data in kibibytes (1024 bytes).
void RemoveReferenceToCell(vtkIdType ptId, vtkIdType cellId)
Remove a reference to a cell in a particular point's link list.
int IsEdge(vtkIdType p1, vtkIdType p2)
Determine whether two points form an edge.
void ReplaceCell(vtkIdType cellId, int npts, const vtkIdType pts[])
Replace the points defining cell "cellId" with a new set of points.
static vtkPolyDataDummyContainter DummyContainer
vtkCellArray * GetPolys()
Get the cell array defining polygons.
void AddCellReference(vtkIdType cellId)
Add references to cell in cell structure.
int GetDataObjectType() override
Return what type of dataset this is.
void ReportReferences(vtkGarbageCollector *) override
vtkSmartPointer< vtkCellArray > Lines
void GetCellEdgeNeighbors(vtkIdType cellId, vtkIdType p1, vtkIdType p2, vtkIdList *cellIds)
Get the neighbors at an edge.
virtual int GetNumberOfPieces()
Get the piece and the number of pieces.
void DeepCopy(vtkDataObject *src) override
Shallow and Deep copy.
void RemoveDeletedCells()
The cells marked by calls to DeleteCell are stored in the Cell Array VTK_EMPTY_CELL,...
vtkIdType GetNumberOfVerts()
Return the number of primitives of a particular type held.
void DeleteCells()
Release data structure that allows random access of the cells.
static vtkPolyData * GetData(vtkInformation *info)
Retrieve an instance of this class from an information object.
int IsPointUsedByCell(vtkIdType ptId, vtkIdType cellId)
Determine whether a point is used by a particular cell.
void ReplaceCellPoint(vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId)
Replace a point in the cell connectivity list with a different point.
vtkSmartPointer< CellMap > Cells
void GetPointCells(vtkIdType ptId, vtkIdList *cellIds) override
Efficient method to obtain cells using a particular point.
void BuildCells()
Create data structure that allows random access of cells.
void DeletePoint(vtkIdType ptId)
Mark a point/cell as deleted from this vtkPolyData.
void CopyCells(vtkPolyData *pd, vtkIdList *idList, vtkIncrementalPointLocator *locator=nullptr)
Copy cells listed in idList from pd, including points, point data, and cell data.
vtkCellArray * GetLines()
Get the cell array defining lines.
void Allocate(vtkIdType numCells=1000, int vtkNotUsed(extSize)=1000)
Method allocates initial storage for vertex, line, polygon, and triangle strip arrays.
void BuildLinks(int initialSize=0)
Create upward links from points to cells that use each point.
vtkIdType GetCellSize(vtkIdType cellId) override
Standard vtkDataSet interface.
void Allocate(vtkPolyData *inPolyData, vtkIdType numCells=1000, int vtkNotUsed(extSize)=1000)
Similar to the method above, this method allocates initial storage for vertex, line,...
void DeleteLinks()
Release the upward links from point to cells that use each point.
vtkMTimeType GetMeshMTime() override
Return the mesh (geometry/topology) modification time.
void ReverseCell(vtkIdType cellId)
Reverse the order of point ids defining the cell.
void SetStrips(vtkCellArray *s)
Set the cell array defining triangle strips.
bool AllocateProportional(vtkPolyData *pd, double ratio)
Preallocate memory for the internal cell arrays such that they are proportional to those in pd by a f...
Hold a reference to a vtkObjectBase instance.
T * Get() const noexcept
Get the contained pointer.
record modification and/or execution time
unsigned char GetCellType() const noexcept
vtkIdType GetCellId() const noexcept
@ VTK_TRIANGLE_STRIP
Definition vtkCellType.h:43
@ VTK_POLY_LINE
Definition vtkCellType.h:41
@ VTK_TRIANGLE
Definition vtkCellType.h:42
@ VTK_POLYGON
Definition vtkCellType.h:44
@ VTK_EMPTY_CELL
Definition vtkCellType.h:37
@ VTK_LINE
Definition vtkCellType.h:40
@ VTK_QUAD
Definition vtkCellType.h:46
@ VTK_VERTEX
Definition vtkCellType.h:38
@ VTK_POLY_VERTEX
Definition vtkCellType.h:39
@ Cells
A cell specified by degrees of freedom held in arrays.
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270
#define VTK_POLY_DATA
Definition vtkType.h:65
#define VTK_SIZEHINT(...)
#define VTK_MARSHALAUTO