VTK  9.3.20240422
vtkPlanesIntersection.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3// SPDX-License-Identifier: BSD-3-Clause
4
49#ifndef vtkPlanesIntersection_h
50#define vtkPlanesIntersection_h
51
52#include "vtkCommonDataModelModule.h" // For export macro
53#include "vtkPlanes.h"
54
55VTK_ABI_NAMESPACE_BEGIN
56class vtkPoints;
58class vtkCell;
59
60class VTKCOMMONDATAMODEL_EXPORT vtkPlanesIntersection : public vtkPlanes
61{
63
64public:
65 void PrintSelf(ostream& os, vtkIndent indent) override;
66
68
76 void SetRegionVertices(double* v, int nvertices);
78 // Retained for backward compatibility
79 int GetNumRegionVertices() { return this->GetNumberOfRegionVertices(); }
80 int GetRegionVertices(double* v, int nvertices);
81
88
97 static int PolygonIntersectsBBox(double bounds[6], vtkPoints* pts);
98
107
108protected:
109 static void ComputeNormal(double* p1, double* p2, double* p3, double normal[3]);
110 static double EvaluatePlaneEquation(double* x, double* p);
111 static void PlaneEquation(double* n, double* x, double* p);
112 static int GoodNormal(double* n);
113 static int Invert3x3(double M[3][3]);
114
117
118private:
119 int IntersectsBoundingBox(vtkPoints* R);
120 int EnclosesBoundingBox(vtkPoints* R);
121 int EvaluateFacePlane(int plane, vtkPoints* R);
122 int IntersectsProjection(vtkPoints* R, int direction);
123
124 void SetPlaneEquations();
125 void ComputeRegionVertices();
126
127 void planesMatrix(int p1, int p2, int p3, double M[3][3]) const;
128 int duplicate(double testv[3]) const;
129 void planesRHS(int p1, int p2, int p3, double r[3]) const;
130 int outsideRegion(double v[3]);
131
132 // plane equations
133 double* Planes;
134
135 // vertices of convex regions enclosed by the planes, also
136 // the ccw hull of that region projected in 3 orthog. directions
137 vtkPointsProjectedHull* RegionPts;
138
140 void operator=(const vtkPlanesIntersection&) = delete;
141};
142VTK_ABI_NAMESPACE_END
143#endif
abstract class to specify cell behavior
Definition vtkCell.h:130
a simple class to control print indentation
Definition vtkIndent.h:108
A vtkPlanesIntersection object is a vtkPlanes object that can compute whether the arbitrary convex re...
static int PolygonIntersectsBBox(double bounds[6], vtkPoints *pts)
A convenience function provided by this class, returns 1 if the polygon defined in pts intersects the...
int GetRegionVertices(double *v, int nvertices)
void SetRegionVertices(double *v, int nvertices)
static void ComputeNormal(double *p1, double *p2, double *p3, double normal[3])
static int GoodNormal(double *n)
static int Invert3x3(double M[3][3])
static double EvaluatePlaneEquation(double *x, double *p)
~vtkPlanesIntersection() override
static vtkPlanesIntersection * Convert3DCell(vtkCell *cell)
Another convenience function provided by this class, returns the vtkPlanesIntersection object represe...
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
void SetRegionVertices(vtkPoints *pts)
It helps if you know the vertices of the convex region.
static vtkPlanesIntersection * New()
static void PlaneEquation(double *n, double *x, double *p)
int IntersectsRegion(vtkPoints *R)
Return 1 if the axis aligned box defined by R intersects the region defined by the planes,...
implicit function for convex set of planes
Definition vtkPlanes.h:151
the convex hull of the orthogonal projection of the vtkPoints in the 3 coordinate directions
represent and manipulate 3D points
Definition vtkPoints.h:139