VTK  9.4.20241121
vtkDataSet.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
143#ifndef vtkDataSet_h
144#define vtkDataSet_h
145
146#include "vtkCommonDataModelModule.h" // For export macro
147#include "vtkDataObject.h"
148#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_3_0
149#include "vtkNew.h" // For vtkNew
150#include "vtkSmartPointer.h" // For vtkSmartPointer
151#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
152
153VTK_ABI_NAMESPACE_BEGIN
154class vtkCell;
155class vtkCellData;
156class vtkCellIterator;
157class vtkCellTypes;
158class vtkGenericCell;
159class vtkIdList;
160class vtkPointData;
161class vtkPoints;
164
165class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALAUTO vtkDataSet : public vtkDataObject
166{
167public:
168 vtkTypeMacro(vtkDataSet, vtkDataObject);
169 void PrintSelf(ostream& os, vtkIndent indent) override;
170
177 virtual void CopyStructure(vtkDataSet* ds) = 0;
178
184 virtual void CopyAttributes(vtkDataSet* ds);
185
191
197
205
210 virtual double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) = 0;
211
218 virtual void GetPoint(vtkIdType id, double x[3]);
219
225
238 virtual vtkCell* GetCell(vtkIdType cellId) = 0;
239 virtual vtkCell* GetCell(int vtkNotUsed(i), int vtkNotUsed(j), int vtkNotUsed(k))
240 {
241 vtkErrorMacro("ijk indices are only valid with structured data!");
242 return nullptr;
243 }
244
246
254 virtual void GetCell(vtkIdType cellId, vtkGenericCell* cell) = 0;
255
267 virtual void GetCellBounds(vtkIdType cellId, double bounds[6]);
268
274 virtual int GetCellType(vtkIdType cellId) = 0;
275
285
295 virtual void GetCellTypes(vtkCellTypes* types);
296
302 virtual void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) = 0;
303
316 virtual void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts,
317 vtkIdList* ptIds) VTK_SIZEHINT(pts, npts);
318
324 virtual void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) = 0;
325
333 virtual void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds);
334
344 int GetCellNumberOfFaces(vtkIdType cellId, unsigned char& cellType, vtkGenericCell* cell);
345
347
354 vtkIdType FindPoint(double x, double y, double z)
355 {
356 double xyz[3];
357 xyz[0] = x;
358 xyz[1] = y;
359 xyz[2] = z;
360 return this->FindPoint(xyz);
361 }
362 virtual vtkIdType FindPoint(double x[3]) = 0;
364
376 virtual vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
377 double pcoords[3], double* weights) = 0;
378
386 virtual vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
387 double tol2, int& subId, double pcoords[3], double* weights) = 0;
388
397 virtual vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2,
398 int& subId, double pcoords[3], double* weights);
399
405
410 vtkCellData* GetCellData() { return this->CellData; }
411
416 vtkPointData* GetPointData() { return this->PointData; }
417
422 virtual void Squeeze();
423
428 virtual void ComputeBounds();
429
436
443 void GetBounds(double bounds[6]);
444
449 double* GetCenter() VTK_SIZEHINT(3);
450
456 void GetCenter(double center[3]);
457
463 double GetLength();
464
470 double GetLength2();
471
476 void Initialize() override;
477
488 virtual void GetScalarRange(double range[2]);
489
499 double* GetScalarRange() VTK_SIZEHINT(2);
500
506 virtual int GetMaxCellSize() = 0;
507
509
516 virtual int GetMaxSpatialDimension();
517 virtual int GetMinSpatialDimension();
519
528 unsigned long GetActualMemorySize() override;
529
533 int GetDataObjectType() override { return VTK_DATA_SET; }
534
536
539 void ShallowCopy(vtkDataObject* src) override;
540 void DeepCopy(vtkDataObject* src) override;
542
544 {
545 DATA_OBJECT_FIELD = 0,
546 POINT_DATA_FIELD = 1,
547 CELL_DATA_FIELD = 2
548 };
549
559
561
566 virtual void GenerateGhostArray(int zeroExt[6]) { this->GenerateGhostArray(zeroExt, false); }
567 virtual void GenerateGhostArray(int zeroExt[6], bool cellOnly);
569
571
577
585
589 vtkIdType GetNumberOfElements(int type) override;
590
601
617 virtual bool HasAnyBlankCells() { return false; }
623 virtual bool HasAnyBlankPoints() { return false; }
624
630
634 VTK_DEPRECATED_IN_9_3_0("This function is deprecated. It has no effect.")
635 void UpdatePointGhostArrayCache() {}
636
641
647
651 VTK_DEPRECATED_IN_9_3_0("This function is deprecated. It has no effect.")
652 void UpdateCellGhostArrayCache() {}
653
664
668 bool SupportsGhostArray(int type) override;
669
670protected:
671 // Constructor with default bounds (0,1, 0,1, 0,1).
673 ~vtkDataSet() override;
674
676
682
687 virtual void ComputeScalarRange();
688
689 vtkCellData* CellData; // Scalars, vectors, etc. associated w/ each cell
690 vtkPointData* PointData; // Scalars, vectors, etc. associated w/ each point
691 vtkCallbackCommand* DataObserver; // Observes changes to cell/point data
692 vtkTimeStamp ComputeTime; // Time at which bounds, center, etc. computed
693 double Bounds[6]; // (xmin,xmax, ymin,ymax, zmin,zmax) geometric bounds
694 double Center[3];
695
696 // Cached scalar range
697 double ScalarRange[2];
698
699 // Time at which scalar range is computed
701
703
707 VTK_DEPRECATED_IN_9_3_0("This member is deprecated. It's no longer used.")
708 vtkUnsignedCharArray* PointGhostArray;
709 VTK_DEPRECATED_IN_9_3_0("This member is deprecated. It's no longer used.")
710 vtkUnsignedCharArray* CellGhostArray;
711 VTK_DEPRECATED_IN_9_3_0("This member is deprecated. It's no longer used.")
712 bool PointGhostArrayCached;
713 VTK_DEPRECATED_IN_9_3_0("This member is deprecated. It's no longer used.")
714 bool CellGhostArrayCached;
716
717private:
718 void InternalDataSetCopy(vtkDataSet* src);
723 static void OnDataModified(
724 vtkObject* source, unsigned long eid, void* clientdata, void* calldata);
725
726 // This should only be used if a vtkDataSet subclass don't define GetPoints()
727 vtkSmartPointer<vtkPoints> TempPoints;
728
729 vtkDataSet(const vtkDataSet&) = delete;
730 void operator=(const vtkDataSet&) = delete;
731};
732
733inline void vtkDataSet::GetPoint(vtkIdType id, double x[3])
734{
735 double* pt = this->GetPoint(id);
736 x[0] = pt[0];
737 x[1] = pt[1];
738 x[2] = pt[2];
739}
740
741VTK_ABI_NAMESPACE_END
742#endif
void GetPoint(int i, int j, int k, double pnt[3])
supports function callbacks
represent and manipulate cell attribute data
Efficient cell iterator for vtkDataSet topologies.
object provides direct access to cells in vtkCellArray and type information
abstract class to specify cell behavior
Definition vtkCell.h:130
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
virtual bool HasAnyBlankPoints()
Returns 1 if there are any blanking points 0 otherwise.
Definition vtkDataSet.h:623
static vtkDataSet * GetData(vtkInformation *info)
Retrieve an instance of this class from an information object.
double * GetBounds()
Return a pointer to the geometry bounding box in the form (xmin,xmax, ymin,ymax, zmin,...
vtkMTimeType GetMTime() override
Datasets are composite objects and need to check each part for MTime THIS METHOD IS THREAD SAFE.
vtkIdType GetNumberOfElements(int type) override
Get the number of elements for a specific attribute type (POINT, CELL, etc.).
virtual vtkIdType GetCellSize(vtkIdType cellId)
Get the size of cell with cellId such that: 0 <= cellId < NumberOfCells.
vtkPointData * GetPointData()
Return a pointer to this dataset's point data.
Definition vtkDataSet.h:416
bool HasAnyGhostCells()
Returns 1 if there are any ghost cells 0 otherwise.
virtual vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)=0
Locate cell based on global coordinate x and tolerance squared.
vtkUnsignedCharArray * AllocateCellGhostArray()
Allocate ghost array for cells.
bool HasAnyGhostPoints()
Returns 1 if there are any ghost points 0 otherwise.
virtual void ComputeBounds()
Compute the data bounding box from data points.
vtkUnsignedCharArray * GetGhostArray(int type) override
Returns the ghost array for the given type (point or cell).
vtkTimeStamp ComputeTime
Definition vtkDataSet.h:692
virtual void CopyAttributes(vtkDataSet *ds)
Copy the attributes associated with the specified dataset to this instance of vtkDataSet.
virtual void GetPointCells(vtkIdType ptId, vtkIdList *cellIds)=0
Topological inquiry to get cells using point.
vtkTimeStamp ScalarRangeComputeTime
Definition vtkDataSet.h:700
static vtkDataSet * GetData(vtkInformationVector *v, int i=0)
Retrieve an instance of this class from an information object.
virtual vtkMTimeType GetMeshMTime()
Abstract method which return the mesh (geometry/topology) modification time.
virtual void GetCellBounds(vtkIdType cellId, double bounds[6])
Get the bounds of the cell with cellId such that: 0 <= cellId < NumberOfCells.
virtual vtkPoints * GetPoints()
If the subclass has (implicit/explicit) points, then return them.
virtual vtkIdType FindCell(double x[3], vtkCell *cell, vtkGenericCell *gencell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)=0
This is a version of the above method that can be used with multithreaded applications.
virtual vtkIdType GetNumberOfPoints()=0
Determine the number of points composing the dataset.
virtual void GetCell(vtkIdType cellId, vtkGenericCell *cell)=0
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
virtual vtkIdType GetNumberOfCells()=0
Determine the number of cells composing the dataset.
vtkFieldData * GetAttributesAsFieldData(int type) override
Returns the attributes of the data object as a vtkFieldData.
virtual bool HasAnyBlankCells()
Returns 1 if there are any blanking cells 0 otherwise.
Definition vtkDataSet.h:617
virtual void GenerateGhostArray(int zeroExt[6])
Normally called by pipeline executives or algorithms only.
Definition vtkDataSet.h:566
vtkPointData * PointData
Definition vtkDataSet.h:690
virtual void ComputeScalarRange()
Compute the range of the scalars and cache it into ScalarRange only if the cache became invalid (Scal...
virtual double * GetPoint(vtkIdType ptId)=0
Get point coordinates with ptId such that: 0 <= ptId < NumberOfPoints.
vtkMTimeType GetGhostCellsTime()
Return the MTime of the ghost cells array.
virtual void GetCellTypes(vtkCellTypes *types)
Get a list of types of cells in a dataset.
vtkNew< vtkGenericCell > GenericCell
Definition vtkDataSet.h:675
virtual void CopyStructure(vtkDataSet *ds)=0
Copy the geometric and topological structure of an object.
virtual vtkCell * GetCell(int vtkNotUsed(i), int vtkNotUsed(j), int vtkNotUsed(k))
Definition vtkDataSet.h:239
vtkCellData * CellData
Definition vtkDataSet.h:689
virtual vtkCell * GetCell(vtkIdType cellId)=0
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
virtual int GetCellType(vtkIdType cellId)=0
Get type of cell with cellId such that: 0 <= cellId < NumberOfCells.
virtual void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)=0
Topological inquiry to get points defining cell.
vtkCallbackCommand * DataObserver
Definition vtkDataSet.h:691
virtual void GetCellPoints(vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts, vtkIdList *ptIds)
Topological inquiry to get points defining cell.
void DeepCopy(vtkDataObject *src) override
Shallow and Deep copy.
vtkUnsignedCharArray * GetCellGhostArray()
Get the array that defines the ghost type of each cell.
int GetCellNumberOfFaces(vtkIdType cellId, unsigned char &cellType, vtkGenericCell *cell)
Get the number of faces of a cell.
vtkUnsignedCharArray * AllocatePointGhostArray()
Allocate ghost array for points.
void SetCellOrderAndRationalWeights(vtkIdType cellId, vtkGenericCell *cell)
vtkIdType FindPoint(double x, double y, double z)
Locate the closest point to the global coordinate x.
Definition vtkDataSet.h:354
virtual vtkCell * FindAndGetCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
Locate the cell that contains a point and return the cell.
virtual void Squeeze()
Reclaim any extra memory used to store data.
~vtkDataSet() override
virtual void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds)
Topological inquiry to get all cells using list of points exclusive of cell specified (e....
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkUnsignedCharArray * GetPointGhostArray()
Gets the array that defines the ghost type of each point.
virtual void GenerateGhostArray(int zeroExt[6], bool cellOnly)
Normally called by pipeline executives or algorithms only.
virtual vtkCellIterator * NewCellIterator()
Return an iterator that traverses the cells in this data set.
int CheckAttributes()
This method checks to see if the cell and point attributes match the geometry.
bool SupportsGhostArray(int type) override
Returns true for POINT or CELL, false otherwise.
vtkCellData * GetCellData()
Return a pointer to this dataset's cell data.
Definition vtkDataSet.h:410
virtual vtkIdType FindPoint(double x[3])=0
Locate the closest point to the global coordinate x.
represent and manipulate fields of data
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:133
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
abstract base class for most VTK objects
Definition vtkObject.h:162
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:139
Hold a reference to a vtkObjectBase instance.
record modification and/or execution time
dynamic, self-adjusting array of unsigned char
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
#define VTK_DEPRECATED_IN_9_3_0(reason)
int vtkIdType
Definition vtkType.h:315
#define VTK_DATA_SET
Definition vtkType.h:73
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270
#define VTK_SIZEHINT(...)
#define VTK_MARSHALAUTO
#define VTK_NEWINSTANCE