VTK  9.3.20240419
vtkFindCellStrategy.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
36 #ifndef vtkFindCellStrategy_h
37 #define vtkFindCellStrategy_h
38 
39 #include "vtkCommonDataModelModule.h" // For export macro
40 #include "vtkObject.h"
41 
42 VTK_ABI_NAMESPACE_BEGIN
43 class vtkCell;
44 class vtkGenericCell;
45 class vtkPointSet;
46 
47 class VTKCOMMONDATAMODEL_EXPORT vtkFindCellStrategy : public vtkObject
48 {
49 public:
51 
55  void PrintSelf(ostream& os, vtkIndent indent) override;
57 
65  virtual int Initialize(vtkPointSet* ps);
66 
75  virtual vtkIdType FindCell(double x[3], vtkCell* cell, vtkGenericCell* gencell, vtkIdType cellId,
76  double tol2, int& subId, double pcoords[3], double* weights) = 0;
77 
94  virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
95  vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) = 0;
96 
100  virtual bool InsideCellBounds(double x[3], vtkIdType cellId) = 0;
101 
110  virtual void CopyParameters(vtkFindCellStrategy* from);
111 
112 protected:
115 
116  // You may ask why this OwnsLocator rigamarole. The reason is that the reference counting garbage
117  // collector gets confused when the (cell/point) locator, point set, and strategy are all mixed
118  // together; resulting in memory leaks etc, So this defines if the locator specified or taken from
119  // another strategy instance or the dataset.
121  // IsACopy is needed to ensure the point-set's locator is up-to-date
122  // otherwise thread-safety issue can arise.
123  bool IsACopy;
124  vtkPointSet* PointSet; // vtkPointSet which this strategy is associated with
125  double Bounds[6]; // bounding box of vtkPointSet
126 
127  vtkTimeStamp InitializeTime; // time at which strategy was initialized
128 
129 private:
130  vtkFindCellStrategy(const vtkFindCellStrategy&) = delete;
131  void operator=(const vtkFindCellStrategy&) = delete;
132 };
133 
134 VTK_ABI_NAMESPACE_END
135 #endif
abstract class to specify cell behavior
Definition: vtkCell.h:130
helper class to manage the vtkPointSet::FindCell() METHOD
virtual bool InsideCellBounds(double x[3], vtkIdType cellId)=0
Quickly test if a point is inside the bounds of a particular cell.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for type information and printing.
virtual void CopyParameters(vtkFindCellStrategy *from)
Copy essential parameters between instances of this class.
virtual vtkIdType FindCell(double x[3], vtkCell *cell, vtkGenericCell *gencell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)=0
Virtual method for finding a cell.
virtual int Initialize(vtkPointSet *ps)
All subclasses of this class must provide an initialize method.
~vtkFindCellStrategy() override
virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2, int &inside)=0
Return the closest point within a specified radius and the cell which is closest to the point x.
provides thread-safe access to cells
a simple class to control print indentation
Definition: vtkIndent.h:108
abstract base class for most VTK objects
Definition: vtkObject.h:162
concrete class for storing a set of points
Definition: vtkPointSet.h:98
record modification and/or execution time
Definition: vtkTimeStamp.h:44
@ radius
Definition: vtkX3D.h:252
int vtkIdType
Definition: vtkType.h:315