VTK  9.0.2
vtkPointSet.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkPointSet.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
51 #ifndef vtkPointSet_h
52 #define vtkPointSet_h
53 
54 #include "vtkCommonDataModelModule.h" // For export macro
55 #include "vtkDataSet.h"
56 
57 #include "vtkPoints.h" // Needed for inline methods
58 
61 
62 class VTKCOMMONDATAMODEL_EXPORT vtkPointSet : public vtkDataSet
63 {
64 public:
66 
69  vtkTypeMacro(vtkPointSet, vtkDataSet);
70  void PrintSelf(ostream& os, vtkIndent indent) override;
71 
73 
82  vtkSetMacro(Editable, bool);
83  vtkGetMacro(Editable, bool);
84  vtkBooleanMacro(Editable, bool);
86 
90  void Initialize() override;
91 
95  void CopyStructure(vtkDataSet* pd) override;
96 
98 
101  vtkIdType GetNumberOfPoints() override;
102  void GetPoint(vtkIdType ptId, double x[3]) override { this->Points->GetPoint(ptId, x); }
103  vtkIdType FindPoint(double x[3]) override;
104  vtkIdType FindPoint(double x, double y, double z) { return this->vtkDataSet::FindPoint(x, y, z); }
105  vtkIdType FindCell(double x[3], vtkCell* cell, vtkIdType cellId, double tol2, int& subId,
106  double pcoords[3], double* weights) override;
107  vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
108  double tol2, int& subId, double pcoords[3], double* weights) override;
110 
117  double* GetPoint(vtkIdType ptId) VTK_SIZEHINT(3) override { return this->Points->GetPoint(ptId); }
118 
123 
125 
130  void BuildLocator() { this->BuildPointLocator(); }
132 
138 
140 
147  vtkGetObjectMacro(PointLocator, vtkAbstractPointLocator);
149 
151 
156  vtkGetObjectMacro(CellLocator, vtkAbstractCellLocator);
158 
162  vtkMTimeType GetMTime() override;
163 
167  void ComputeBounds() override;
168 
172  void Squeeze() override;
173 
175 
178  virtual void SetPoints(vtkPoints*);
179  vtkGetObjectMacro(Points, vtkPoints);
181 
190  unsigned long GetActualMemorySize() override;
191 
193 
196  void ShallowCopy(vtkDataObject* src) override;
197  void DeepCopy(vtkDataObject* src) override;
199 
201 
204  void Register(vtkObjectBase* o) override;
205  void UnRegister(vtkObjectBase* o) override;
207 
209 
213  static vtkPointSet* GetData(vtkInformationVector* v, int i = 0);
215 
216 protected:
218  ~vtkPointSet() override;
219 
220  bool Editable;
224 
226 
227 private:
228  void Cleanup();
229 
230  vtkPointSet(const vtkPointSet&) = delete;
231  void operator=(const vtkPointSet&) = delete;
232 };
233 
235 {
236  if (this->Points)
237  {
238  return this->Points->GetNumberOfPoints();
239  }
240  else
241  {
242  return 0;
243  }
244 }
245 
246 #endif
an abstract base class for locators which find cells
abstract class to quickly locate points in 3-space
Efficient cell iterator for vtkDataSet topologies.
abstract class to specify cell behavior
Definition: vtkCell.h:57
general representation of visualization data
Definition: vtkDataObject.h:60
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
virtual vtkIdType GetNumberOfPoints()=0
Determine the number of points composing the dataset.
vtkIdType FindPoint(double x, double y, double z)
Locate the closest point to the global coordinate x.
Definition: vtkDataSet.h:193
Detect and break reference loops.
provides thread-safe access to cells
a simple class to control print indentation
Definition: vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
abstract base class for most VTK objects
Definition: vtkObjectBase.h:64
abstract class for specifying dataset behavior
Definition: vtkPointSet.h:63
void CopyStructure(vtkDataSet *pd) override
Copy the geometric structure of an input point set object.
void Initialize() override
Reset to an empty state and free any memory.
vtkIdType FindPoint(double x, double y, double z)
Definition: vtkPointSet.h:104
virtual void SetCellLocator(vtkAbstractCellLocator *)
Set / get an instance of vtkAbstractCellLocator which may be used when a vtkCellLocatorStrategy is us...
void BuildPointLocator()
Build the internal point locator .
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual void SetPoints(vtkPoints *)
Specify point array to define point coordinates.
~vtkPointSet() override
virtual void SetPointLocator(vtkAbstractPointLocator *)
Set / get an instance of vtkAbstractPointLocator which is used to support the FindPoint() and FindCel...
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
void Squeeze() override
Reclaim any unused memory.
static vtkPointSet * GetData(vtkInformationVector *v, int i=0)
vtkIdType FindPoint(double x[3]) override
vtkIdType FindCell(double x[3], vtkCell *cell, vtkGenericCell *gencell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
This is a version of the above method that can be used with multithreaded applications.
vtkCellIterator * NewCellIterator() override
Return an iterator that traverses the cells in this data set.
vtkAbstractPointLocator * PointLocator
Definition: vtkPointSet.h:222
void UnRegister(vtkObjectBase *o) override
Decrease the reference count (release by another object).
void Register(vtkObjectBase *o) override
Overwritten to handle the data/locator loop.
void BuildLocator()
Definition: vtkPointSet.h:130
vtkIdType FindCell(double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
Locate cell based on global coordinate x and tolerance squared.
vtkAbstractCellLocator * CellLocator
Definition: vtkPointSet.h:223
void BuildCellLocator()
Build the cell locator.
void ComputeBounds() override
Compute the (X, Y, Z) bounds of the data.
vtkPoints * Points
Definition: vtkPointSet.h:221
vtkIdType GetNumberOfPoints() override
See vtkDataSet for additional information.
Definition: vtkPointSet.h:234
double * GetPoint(vtkIdType ptId) override
See vtkDataSet for additional information.
Definition: vtkPointSet.h:117
vtkMTimeType GetMTime() override
Get MTime which also considers its vtkPoints MTime.
void DeepCopy(vtkDataObject *src) override
void GetPoint(vtkIdType ptId, double x[3]) override
Copy point coordinates into user provided array x[3] for specified point id.
Definition: vtkPointSet.h:102
void ReportReferences(vtkGarbageCollector *) override
static vtkPointSet * GetData(vtkInformation *info)
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).
represent and manipulate 3D points
Definition: vtkPoints.h:34
vtkIdType GetNumberOfPoints()
Return number of points in array.
Definition: vtkPoints.h:125
@ info
Definition: vtkX3D.h:382
int vtkIdType
Definition: vtkType.h:338
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:293
#define VTK_SIZEHINT(...)