VTK  9.6.20260602
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
142
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_6_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;
157class vtkCellTypes;
158class vtkGenericCell;
159class vtkIdList;
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
287
296 virtual void GetDistinctCellTypes(vtkCellTypes* types);
297 VTK_DEPRECATED_IN_9_6_0("Use GetDistinctCellTypes(vtkCellTypes* types) instead.")
298 virtual void GetCellTypes(vtkCellTypes* types) { this->GetDistinctCellTypes(types); }
300
306 virtual void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) = 0;
307
320 virtual void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts,
321 vtkIdList* ptIds) VTK_SIZEHINT(pts, npts);
322
328 virtual void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) = 0;
329
337 virtual void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds);
338
347 virtual int GetCellNumberOfFaces(vtkIdType cellId, unsigned char& cellType, vtkGenericCell* cell);
348
350
357 vtkIdType FindPoint(double x, double y, double z)
358 {
359 double xyz[3];
360 xyz[0] = x;
361 xyz[1] = y;
362 xyz[2] = z;
363 return this->FindPoint(xyz);
364 }
365 virtual vtkIdType FindPoint(double x[3]) = 0;
367
379 virtual vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
380 double pcoords[3], double* weights);
381
389 virtual vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
390 double tol2, int& subId, double pcoords[3], double* weights) = 0;
391
400 virtual vtkCell* FindAndGetCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2,
401 int& subId, double pcoords[3], double* weights);
402
408
413 vtkCellData* GetCellData() { return this->CellData; }
414
420
425 virtual void Squeeze();
426
431 virtual void ComputeBounds();
432
439
446 void GetBounds(double bounds[6]);
447
453
459 void GetCenter(double center[3]);
460
466 double GetLength();
467
473 double GetLength2();
474
490 double GetSampledMaxCellLength2(vtkIdType numSamples = 100);
491
496 void Initialize() override;
497
508 virtual void GetScalarRange(double range[2]);
509
520
526 virtual int GetMaxCellSize() = 0;
527
529
539
548 unsigned long GetActualMemorySize() override;
549
553 int GetDataObjectType() VTK_FUTURE_CONST override { return VTK_DATA_SET; }
554
556
559 void ShallowCopy(vtkDataObject* src) override;
560 void DeepCopy(vtkDataObject* src) override;
562
569
579
581
586 virtual void GenerateGhostArray(int zeroExt[6]) { this->GenerateGhostArray(zeroExt, false); }
587 virtual void GenerateGhostArray(int zeroExt[6], bool cellOnly);
589
591
597
605
609 vtkIdType GetNumberOfElements(int type) override;
610
621
637 virtual bool HasAnyBlankCells() { return false; }
643 virtual bool HasAnyBlankPoints() { return false; }
644
650
655
661
672
676 bool SupportsGhostArray(int type) override;
677
678protected:
679 // Constructor with default bounds (0,1, 0,1, 0,1).
681 ~vtkDataSet() override;
682
684
690
695 virtual void ComputeScalarRange();
696
697 vtkCellData* CellData; // Scalars, vectors, etc. associated w/ each cell
698 vtkPointData* PointData; // Scalars, vectors, etc. associated w/ each point
699 vtkCallbackCommand* DataObserver; // Observes changes to cell/point data
700 vtkTimeStamp ComputeTime; // Time at which bounds, center, etc. computed
701 double Bounds[6]; // (xmin,xmax, ymin,ymax, zmin,zmax) geometric bounds
702 double Center[3];
703
704 // Cached scalar range
705 double ScalarRange[2];
706
707 // Time at which scalar range is computed
709
710private:
711 void InternalDataSetCopy(vtkDataSet* src);
712
713 // This should only be used if a vtkDataSet subclass don't define GetPoints()
715
716 vtkDataSet(const vtkDataSet&) = delete;
717 void operator=(const vtkDataSet&) = delete;
718};
719
720inline void vtkDataSet::GetPoint(vtkIdType id, double x[3])
721{
722 double* pt = this->GetPoint(id);
723 x[0] = pt[0];
724 x[1] = pt[1];
725 x[2] = pt[2];
726}
727
728VTK_ABI_NAMESPACE_END
729#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
virtual vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
Locate cell based on global coordinate x and tolerance squared.
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:643
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:419
bool HasAnyGhostCells()
Returns 1 if there are any ghost cells 0 otherwise.
double Center[3]
Definition vtkDataSet.h:702
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:700
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:708
static vtkDataSet * GetData(vtkInformationVector *v, int i=0)
Retrieve an instance of this class from an information object.
unsigned long GetActualMemorySize() override
Return the actual size of the data in kibibytes (1024 bytes).
int GetDataObjectType() VTK_FUTURE_CONST override
Return the type of data object.
Definition vtkDataSet.h:553
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.
double GetLength()
Return the length of the diagonal of the bounding box.
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.
double GetLength2()
Return the squared length of the diagonal of the bounding box.
virtual vtkIdType GetNumberOfPoints()=0
Determine the number of points composing the dataset.
@ DATA_OBJECT_FIELD
Definition vtkDataSet.h:565
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:637
virtual void GenerateGhostArray(int zeroExt[6])
Normally called by pipeline executives or algorithms only.
Definition vtkDataSet.h:586
virtual int GetMaxSpatialDimension()
Get the maximum/minimum spatial dimensionality of the data which is the maximum/minimum dimension of ...
vtkPointData * PointData
Definition vtkDataSet.h:698
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.
Definition vtkDataSet.h:298
virtual void GetScalarRange(double range[2])
Convenience method to get the range of the first component (and only the first component) of any scal...
double * GetCenter()
Get the center of the bounding box.
vtkNew< vtkGenericCell > GenericCell
Definition vtkDataSet.h:683
virtual void CopyStructure(vtkDataSet *ds)=0
Copy the geometric and topological structure of an object.
void Initialize() override
Restore data object to initial state.
double Bounds[6]
Definition vtkDataSet.h:701
virtual void GetDistinctCellTypes(vtkCellTypes *types)
Get a list of types of cells in a dataset.
double GetSampledMaxCellLength2(vtkIdType numSamples=100)
Samples max squared cell length using strided indexing.
vtkCellData * CellData
Definition vtkDataSet.h:697
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 int GetMinSpatialDimension()
Get the maximum/minimum spatial dimensionality of the data which is the maximum/minimum dimension of ...
virtual void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)=0
Topological inquiry to get points defining cell.
vtkCallbackCommand * DataObserver
Definition vtkDataSet.h:699
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.
virtual 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:357
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.
double ScalarRange[2]
Definition vtkDataSet.h:705
~vtkDataSet() override
virtual int GetMaxCellSize()=0
Convenience method returns largest cell size in dataset.
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 vtkCell * GetCell(int i, int j, int k)
Definition vtkDataSet.h:239
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:413
virtual vtkIdType FindPoint(double x[3])=0
Locate the closest point to the global coordinate x.
Represents and manipulates a collection of data arrays.
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:135
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:168
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:140
Hold a reference to a vtkObjectBase instance.
record modification and/or execution time
dynamic, self-adjusting array of unsigned char
#define VTK_DEPRECATED_IN_9_6_0(reason)
int vtkIdType
Definition vtkType.h:363
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:318
@ VTK_DATA_SET
Definition vtkType.h:117
#define VTK_SIZEHINT(...)
#define VTK_MARSHALAUTO
#define VTK_NEWINSTANCE