VTK  9.4.20251007
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
55
56#ifndef vtkPolyData_h
57#define vtkPolyData_h
58
59#include "vtkCommonDataModelModule.h" // For export macro
60#include "vtkPointSet.h"
61#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
62
63#include "vtkCellArray.h" // Needed for inline methods
64#include "vtkCellLinks.h" // Needed for inline methods
65#include "vtkPolyDataInternals.h" // Needed for inline methods
66
67VTK_ABI_NAMESPACE_BEGIN
68struct vtkPolyDataDummyContainter;
70
71class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALAUTO vtkPolyData : public vtkPointSet
72{
73public:
74 static vtkPolyData* New();
76
77 vtkTypeMacro(vtkPolyData, vtkPointSet);
78 void PrintSelf(ostream& os, vtkIndent indent) override;
79
83 int GetDataObjectType() override { return VTK_POLY_DATA; }
84
88 void CopyStructure(vtkDataSet* ds) override;
89
91
94 vtkIdType GetNumberOfCells() override;
96 vtkCell* GetCell(vtkIdType cellId) override;
97 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
98 int GetCellType(vtkIdType cellId) override;
99 vtkIdType GetCellSize(vtkIdType cellId) override;
100 void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
101 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override;
103
111 void CopyCells(vtkPolyData* pd, vtkIdList* idList, vtkIncrementalPointLocator* locator = nullptr);
112
116 void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override;
117
122 void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override;
123
143
149 void GetCellsBounds(double bounds[6]);
150
157 void Squeeze() override;
158
162 int GetMaxCellSize() override;
163
169
176
181
187
192
198
203
209
214
221
223
226 vtkIdType GetNumberOfVerts() { return (this->Verts ? this->Verts->GetNumberOfCells() : 0); }
227 vtkIdType GetNumberOfLines() { return (this->Lines ? this->Lines->GetNumberOfCells() : 0); }
228 vtkIdType GetNumberOfPolys() { return (this->Polys ? this->Polys->GetNumberOfCells() : 0); }
229 vtkIdType GetNumberOfStrips() { return (this->Strips ? this->Strips->GetNumberOfCells() : 0); }
231
241 bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize);
242
252 bool AllocateEstimate(vtkIdType numVerts, vtkIdType maxVertSize, vtkIdType numLines,
253 vtkIdType maxLineSize, vtkIdType numPolys, vtkIdType maxPolySize, vtkIdType numStrips,
254 vtkIdType maxStripSize);
255
265 bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize);
266
277 bool AllocateExact(vtkIdType numVerts, vtkIdType vertConnSize, vtkIdType numLines,
278 vtkIdType lineConnSize, vtkIdType numPolys, vtkIdType polyConnSize, vtkIdType numStrips,
279 vtkIdType stripConnSize);
280
290
300 bool AllocateProportional(vtkPolyData* pd, double ratio);
301
308 void Allocate(vtkIdType numCells = 1000, int vtkNotUsed(extSize) = 1000)
309 {
310 this->AllocateExact(numCells, numCells);
311 }
312
323 void Allocate(vtkPolyData* inPolyData, vtkIdType numCells = 1000, int vtkNotUsed(extSize) = 1000)
324 {
326 inPolyData, static_cast<double>(numCells) / inPolyData->GetNumberOfCells());
327 }
328
336 vtkIdType InsertNextCell(int type, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
337
346
351 void Reset();
352
361
365 bool NeedToBuildCells() { return this->Cells == nullptr; }
366
373 void BuildLinks(int initialSize = 0);
374
376
379 vtkSetSmartPointerMacro(Links, vtkAbstractCellLinks);
380 vtkGetSmartPointerMacro(Links, vtkAbstractCellLinks);
382
389
394
396
400 void GetPointCells(vtkIdType ptId, vtkIdType& ncells, vtkIdType*& cells)
401 VTK_SIZEHINT(cells, ncells);
403
410
422 unsigned char GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts)
423 VTK_SIZEHINT(pts, npts);
424
439 void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds)
440 VTK_SIZEHINT(pts, npts) override;
441
446 int IsTriangle(int v1, int v2, int v3);
447
456
461 int IsPointUsedByCell(vtkIdType ptId, vtkIdType cellId);
462
471 void ReplaceCell(vtkIdType cellId, vtkIdList* ids);
472 void ReplaceCell(vtkIdType cellId, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
474
476
485 void ReplaceCellPoint(vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId);
486 void ReplaceCellPoint(
487 vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId, vtkIdList* cellPointIds);
489
494 void ReverseCell(vtkIdType cellId);
495
497
501 void DeletePoint(vtkIdType ptId);
502 void DeleteCell(vtkIdType cellId);
504
514
516
525 vtkIdType InsertNextLinkedPoint(double x[3], int numLinks);
527
534 vtkIdType InsertNextLinkedCell(int type, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
535
545 void ReplaceLinkedCell(vtkIdType cellId, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
546
554 void RemoveCellReference(vtkIdType cellId);
555
563 void AddCellReference(vtkIdType cellId);
564
573
582
588 void ResizeCellList(vtkIdType ptId, int size);
589
593 void Initialize() override;
594
596
599 virtual int GetPiece();
600 virtual int GetNumberOfPieces();
602
606 virtual int GetGhostLevel();
607
616 unsigned long GetActualMemorySize() override;
617
619
622 void ShallowCopy(vtkDataObject* src) override;
623 void DeepCopy(vtkDataObject* src) override;
625
633
635
641
660 enum
661 {
669 };
670
672 int GetScalarFieldCriticalIndex(vtkIdType pointId, int fieldId);
673 int GetScalarFieldCriticalIndex(vtkIdType pointId, const char* fieldName);
674
683
688
698 unsigned char GetCell(vtkIdType cellId, const vtkIdType*& pts);
699
700protected:
702 ~vtkPolyData() override;
703
705
708
710
711 // points inherited
712 // point data (i.e., scalars, vectors, normals, tcoords) inherited
717
718 // supporting structures for more complex topological operations
719 // built only when necessary
722
724
725 // dummy static member below used as a trick to simplify traversal
726 static vtkPolyDataDummyContainter DummyContainer;
727
728 // Take into account only points that belong to at least one cell.
729 double CellsBounds[6];
730
732
733private:
734 void Cleanup();
735
736 vtkPolyData(const vtkPolyData&) = delete;
737 void operator=(const vtkPolyData&) = delete;
738};
739
740//------------------------------------------------------------------------------
742{
743 return (this->GetNumberOfVerts() + this->GetNumberOfLines() + this->GetNumberOfPolys() +
744 this->GetNumberOfStrips());
745}
746
747//------------------------------------------------------------------------------
749{
750 if (!this->Cells)
751 {
752 this->BuildCells();
753 }
754 return static_cast<int>(this->Cells->GetTag(cellId).GetCellType());
755}
756
757//------------------------------------------------------------------------------
759{
760 if (!this->Cells)
761 {
762 this->BuildCells();
763 }
764 switch (this->GetCellType(cellId))
765 {
766 case VTK_EMPTY_CELL:
767 return 0;
768 case VTK_VERTEX:
769 return 1;
770 case VTK_LINE:
771 return 2;
772 case VTK_TRIANGLE:
773 return 3;
774 case VTK_QUAD:
775 return 4;
776 case VTK_POLY_VERTEX:
777 return this->Verts ? this->Verts->GetCellSize(this->GetCellIdRelativeToCellArray(cellId)) : 0;
778 case VTK_POLY_LINE:
779 return this->Lines ? this->Lines->GetCellSize(this->GetCellIdRelativeToCellArray(cellId)) : 0;
780 case VTK_POLYGON:
781 return this->Polys ? this->Polys->GetCellSize(this->GetCellIdRelativeToCellArray(cellId)) : 0;
783 return this->Strips ? this->Strips->GetCellSize(this->GetCellIdRelativeToCellArray(cellId))
784 : 0;
785 }
786 vtkWarningMacro(<< "Cell type not supported.");
787 return 0;
788}
789
790//------------------------------------------------------------------------------
792{
793 vtkIdType npts;
794 const vtkIdType* pts;
795
796 this->GetCellPoints(cellId, npts, pts);
797 for (vtkIdType i = 0; i < npts; i++)
798 {
799 if (pts[i] == ptId)
800 {
801 return 1;
802 }
803 }
804
805 return 0;
806}
807
808//------------------------------------------------------------------------------
810{
811 static_cast<vtkCellLinks*>(this->Links.Get())->DeletePoint(ptId);
812}
813
814//------------------------------------------------------------------------------
816{
817 this->Cells->GetTag(cellId).MarkDeleted();
818}
819
820//------------------------------------------------------------------------------
822{
823 const vtkIdType* pts;
824 vtkIdType npts;
825
826 this->GetCellPoints(cellId, npts, pts);
827 auto links = static_cast<vtkCellLinks*>(this->Links.Get());
828 for (vtkIdType i = 0; i < npts; i++)
829 {
830 links->RemoveCellReference(cellId, pts[i]);
831 }
832}
833
834//------------------------------------------------------------------------------
836{
837 const vtkIdType* pts;
838 vtkIdType npts;
839
840 this->GetCellPoints(cellId, npts, pts);
841 auto links = static_cast<vtkCellLinks*>(this->Links.Get());
842 for (vtkIdType i = 0; i < npts; i++)
843 {
844 links->AddCellReference(cellId, pts[i]);
845 }
846}
847
848//------------------------------------------------------------------------------
849inline void vtkPolyData::ResizeCellList(vtkIdType ptId, int size)
850{
851 static_cast<vtkCellLinks*>(this->Links.Get())->ResizeCellList(ptId, size);
852}
853
854//------------------------------------------------------------------------------
856{
857 switch (tag.GetTarget())
858 {
860 return this->Verts;
862 return this->Lines;
864 return this->Polys;
866 return this->Strips;
867 }
868 return nullptr; // unreachable
869}
870
871//------------------------------------------------------------------------------
872inline void vtkPolyData::ReplaceCellPoint(vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId)
873{
875 this->ReplaceCellPoint(cellId, oldPtId, newPtId, ids);
876}
877
878//------------------------------------------------------------------------------
880 vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId, vtkIdList* cellPointIds)
881{
882 if (!this->Cells)
883 {
884 this->BuildCells();
885 }
886 vtkIdType npts;
887 const vtkIdType* pts;
888 this->GetCellPoints(cellId, npts, pts, cellPointIds);
889 for (vtkIdType i = 0; i < npts; i++)
890 {
891 if (pts[i] == oldPtId)
892 {
893 const TaggedCellId tag = this->Cells->GetTag(cellId);
894 vtkCellArray* cells = this->GetCellArrayInternal(tag);
895 cells->ReplaceCellPointAtId(tag.GetCellId(), i, newPtId);
896 break;
897 }
898 }
899}
900
901//------------------------------------------------------------------------------
902inline unsigned char vtkPolyData::GetCellPoints(
903 vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts)
904{
905 if (!this->Cells)
906 {
907 this->BuildCells();
908 }
909
910 const TaggedCellId tag = this->Cells->GetTag(cellId);
911 if (tag.IsDeleted())
912 {
913 npts = 0;
914 pts = nullptr;
915 return VTK_EMPTY_CELL;
916 }
917
918 vtkCellArray* cells = this->GetCellArrayInternal(tag);
919 cells->GetCellAtId(tag.GetCellId(), npts, pts);
920 return tag.GetCellType();
921}
922
923//------------------------------------------------------------------------------
925 vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds)
926{
927 if (!this->Cells)
928 {
929 this->BuildCells();
930 }
931
932 const TaggedCellId tag = this->Cells->GetTag(cellId);
933 if (tag.IsDeleted())
934 {
935 npts = 0;
936 pts = nullptr;
937 }
938 else
939 {
940 vtkCellArray* cells = this->GetCellArrayInternal(tag);
941 cells->GetCellAtId(tag.GetCellId(), npts, pts, ptIds);
942 }
943}
944
945VTK_ABI_NAMESPACE_END
946#endif
object to represent cell connectivity
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:51
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:24
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Allocate and hold a VTK object.
Definition vtkNew.h:58
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.
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.
double CellsBounds[6]
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.
void CopyStructure(vtkDataSet *ds) override
Copy the geometric and topological structure of an input poly data object.
vtkNew< vtkIdList > LegacyBuffer
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.
vtkPolyData_detail::TaggedCellId TaggedCellId
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.
@ ERR_NON_MANIFOLD_STAR
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...
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 spatial dimensionality of the data which is the maximum dimension of all cells.
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.
Definition vtkPolyData.h:83
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.
vtkPolyData_detail::CellMap CellMap
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.
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
#define vtkDataArray
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