|
VTK
9.0.1
|
a 3D cell defined by a set of polygonal faces More...
#include <vtkPolyhedron.h>
Public Types | |
| typedef vtkCell3D | Superclass |
Public Types inherited from vtkCell3D | |
| typedef vtkCell | Superclass |
Public Types inherited from vtkCell | |
| typedef vtkObject | Superclass |
Public Member Functions | |
| virtual vtkTypeBool | IsA (const char *type) |
| Return 1 if this class is the same type of (or a subclass of) the named class. More... | |
| vtkPolyhedron * | NewInstance () const |
| void | PrintSelf (ostream &os, vtkIndent indent) override |
| Methods invoked by print to print information about the object including superclasses. More... | |
| void | GetEdgePoints (vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override |
| See vtkCell3D API for description of these methods. More... | |
| void | GetEdgePoints (int vtkNotUsed(edgeId), int *&vtkNotUsed(pts)) override |
| vtkIdType | GetFacePoints (vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(pts)) override |
| void | GetFacePoints (int vtkNotUsed(faceId), int *&vtkNotUsed(pts)) override |
| void | GetEdgeToAdjacentFaces (vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override |
| vtkIdType | GetFaceToAdjacentFaces (vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(faceIds)) override |
| vtkIdType | GetPointToIncidentEdges (vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(edgeIds)) override |
| vtkIdType | GetPointToIncidentFaces (vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(faceIds)) override |
| vtkIdType | GetPointToOneRingPoints (vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(pts)) override |
| bool | GetCentroid (double vtkNotUsed(centroid)[3]) const override |
| double * | GetParametricCoords () override |
| See vtkCell3D API for description of this method. More... | |
| int | GetCellType () override |
| See the vtkCell API for descriptions of these methods. More... | |
| int | RequiresInitialization () override |
| This cell requires that it be initialized prior to access. More... | |
| void | Initialize () override |
| int | GetNumberOfEdges () override |
| A polyhedron is represented internally by a set of polygonal faces. More... | |
| vtkCell * | GetEdge (int) override |
| Return the edge cell from the edgeId of the cell. More... | |
| int | GetNumberOfFaces () override |
| Return the number of faces in the cell. More... | |
| vtkCell * | GetFace (int faceId) override |
| Return the face cell from the faceId of the cell. More... | |
| void | Contour (double value, vtkDataArray *scalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override |
| Satisfy the vtkCell API. More... | |
| void | Clip (double value, vtkDataArray *scalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override |
| Satisfy the vtkCell API. More... | |
| int | EvaluatePosition (const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override |
| Satisfy the vtkCell API. More... | |
| void | EvaluateLocation (int &subId, const double pcoords[3], double x[3], double *weights) override |
| The inverse of EvaluatePosition. More... | |
| int | IntersectWithLine (const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override |
| Intersect the line (p1,p2) with a given tolerance tol to determine a point of intersection x[3] with parametric coordinate t along the line. More... | |
| int | Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts) override |
| Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh. More... | |
| void | Derivatives (int subId, const double pcoords[3], const double *values, int dim, double *derivs) override |
| Computes derivatives at the point specified by the parameter coordinate. More... | |
| int | CellBoundary (int subId, const double pcoords[3], vtkIdList *pts) override |
| Find the boundary face closest to the point defined by the pcoords[3] and subId of the cell (subId can be ignored). More... | |
| int | GetParametricCenter (double pcoords[3]) override |
| Return the center of the cell in parametric coordinates. More... | |
| int | IsPrimaryCell () override |
| A polyhedron is a full-fledged primary cell. More... | |
| void | InterpolateFunctions (const double x[3], double *sf) override |
| Compute the interpolation functions/derivatives (aka shape functions/derivatives). More... | |
| void | InterpolateDerivs (const double x[3], double *derivs) override |
| int | RequiresExplicitFaceRepresentation () override |
| Methods supporting the definition of faces. More... | |
| void | SetFaces (vtkIdType *faces) override |
| vtkIdType * | GetFaces () override |
| int | IsInside (const double x[3], double tolerance) |
| A method particular to vtkPolyhedron. More... | |
| bool | IsConvex () |
| Determine whether or not a polyhedron is convex. More... | |
| vtkPolyData * | GetPolyData () |
| Construct polydata if no one exist, then return this->PolyData. More... | |
Public Member Functions inherited from vtkCell3D | |
| vtkCell3D * | NewInstance () const |
| virtual void | GetEdgePoints (vtkIdType edgeId, const vtkIdType *&pts)=0 |
| Get the pair of vertices that define an edge. More... | |
| virtual void | GetEdgePoints (int edgeId, int *&pts)=0 |
| virtual vtkIdType | GetFacePoints (vtkIdType faceId, const vtkIdType *&pts)=0 |
| Get the list of vertices that define a face. More... | |
| virtual void | GetFacePoints (int faceId, int *&pts)=0 |
| virtual void | GetEdgeToAdjacentFaces (vtkIdType edgeId, const vtkIdType *&faceIds)=0 |
| Get the ids of the two adjacent faces to edge of id edgeId. More... | |
| virtual vtkIdType | GetFaceToAdjacentFaces (vtkIdType faceId, const vtkIdType *&faceIds)=0 |
| Get the ids of the adjacent faces to face of id faceId. More... | |
| virtual vtkIdType | GetPointToIncidentEdges (vtkIdType pointId, const vtkIdType *&edgeIds)=0 |
| Get the ids of the incident edges to point of id pointId. More... | |
| virtual vtkIdType | GetPointToIncidentFaces (vtkIdType pointId, const vtkIdType *&faceIds)=0 |
| Get the ids of the incident faces point of id pointId. More... | |
| virtual vtkIdType | GetPointToOneRingPoints (vtkIdType pointId, const vtkIdType *&pts)=0 |
| Get the ids of a one-ring surrounding point of id pointId. More... | |
| virtual bool | IsInsideOut () |
| Returns true if the normals of the vtkCell3D point inside the cell. More... | |
| virtual bool | GetCentroid (double centroid[3]) const =0 |
| Computes the centroid of the cell. More... | |
| int | GetCellDimension () override |
| The topological dimension of the cell. More... | |
| virtual void | SetMergeTolerance (double) |
| Set the tolerance for merging clip intersection points that are near the vertices of cells. More... | |
| virtual double | GetMergeTolerance () |
Public Member Functions inherited from vtkCell | |
| vtkCell * | NewInstance () const |
| void | Initialize (int npts, const vtkIdType *pts, vtkPoints *p) |
| Initialize cell from outside with point ids and point coordinates specified. More... | |
| void | Initialize (int npts, vtkPoints *p) |
| Initialize the cell with point coordinates specified. More... | |
| virtual void | ShallowCopy (vtkCell *c) |
| Copy this cell by reference counting the internal data structures. More... | |
| virtual void | DeepCopy (vtkCell *c) |
| Copy this cell by completely copying internal data structures. More... | |
| virtual int | IsLinear () |
| Non-linear cells require special treatment beyond the usual cell type and connectivity list information. More... | |
| virtual int | IsExplicitCell () |
| Explicit cells require additional representational information beyond the usual cell type and connectivity list information. More... | |
| virtual void | SetFaces (vtkIdType *vtkNotUsed(faces)) |
| vtkPoints * | GetPoints () |
| Get the point coordinates for the cell. More... | |
| vtkIdType | GetNumberOfPoints () |
| Return the number of points in the cell. More... | |
| vtkIdList * | GetPointIds () |
| Return the list of point ids defining the cell. More... | |
| vtkIdType | GetPointId (int ptId) |
| For cell point i, return the actual point id. More... | |
| void | GetBounds (double bounds[6]) |
| Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More... | |
| double * | GetBounds () |
| Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More... | |
| double | GetLength2 () |
| Compute Length squared of cell (i.e., bounding box diagonal squared). More... | |
| virtual double | GetParametricDistance (const double pcoords[3]) |
| Return the distance of the parametric coordinate provided to the cell. More... | |
| virtual void | InterpolateFunctions (const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight)) |
| Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this level. More... | |
| virtual void | InterpolateDerivs (const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs)) |
Public Member Functions inherited from vtkObject | |
| vtkBaseTypeMacro (vtkObject, vtkObjectBase) | |
| virtual void | DebugOn () |
| Turn debugging output on. More... | |
| virtual void | DebugOff () |
| Turn debugging output off. More... | |
| bool | GetDebug () |
| Get the value of the debug flag. More... | |
| void | SetDebug (bool debugFlag) |
| Set the value of the debug flag. More... | |
| virtual void | Modified () |
| Update the modification time for this object. More... | |
| virtual vtkMTimeType | GetMTime () |
| Return this object's modified time. More... | |
| unsigned long | AddObserver (unsigned long event, vtkCommand *, float priority=0.0f) |
| Allow people to add/remove/invoke observers (callbacks) to any VTK object. More... | |
| unsigned long | AddObserver (const char *event, vtkCommand *, float priority=0.0f) |
| vtkCommand * | GetCommand (unsigned long tag) |
| void | RemoveObserver (vtkCommand *) |
| void | RemoveObservers (unsigned long event, vtkCommand *) |
| void | RemoveObservers (const char *event, vtkCommand *) |
| vtkTypeBool | HasObserver (unsigned long event, vtkCommand *) |
| vtkTypeBool | HasObserver (const char *event, vtkCommand *) |
| void | RemoveObserver (unsigned long tag) |
| void | RemoveObservers (unsigned long event) |
| void | RemoveObservers (const char *event) |
| void | RemoveAllObservers () |
| vtkTypeBool | HasObserver (unsigned long event) |
| vtkTypeBool | HasObserver (const char *event) |
| template<class U , class T > | |
| unsigned long | AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f) |
| Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More... | |
| template<class U , class T > | |
| unsigned long | AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f) |
| template<class U , class T > | |
| unsigned long | AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f) |
| Allow user to set the AbortFlagOn() with the return value of the callback method. More... | |
| int | InvokeEvent (unsigned long event, void *callData) |
| This method invokes an event and return whether the event was aborted or not. More... | |
| int | InvokeEvent (const char *event, void *callData) |
| int | InvokeEvent (unsigned long event) |
| int | InvokeEvent (const char *event) |
Public Member Functions inherited from vtkObjectBase | |
| const char * | GetClassName () const |
| Return the class name as a string. More... | |
| virtual vtkIdType | GetNumberOfGenerationsFromBase (const char *name) |
| Given a the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class). More... | |
| virtual void | Delete () |
| Delete a VTK object. More... | |
| virtual void | FastDelete () |
| Delete a reference to this object. More... | |
| void | InitializeObjectBase () |
| void | Print (ostream &os) |
| Print an object to an ostream. More... | |
| virtual void | PrintHeader (ostream &os, vtkIndent indent) |
| virtual void | PrintTrailer (ostream &os, vtkIndent indent) |
| virtual void | Register (vtkObjectBase *o) |
| Increase the reference count (mark as used by another object). More... | |
| virtual void | UnRegister (vtkObjectBase *o) |
| Decrease the reference count (release by another object). More... | |
| int | GetReferenceCount () |
| Return the current reference count of this object. More... | |
| void | SetReferenceCount (int) |
| Sets the reference count. More... | |
| void | PrintRevisions (ostream &) |
| Legacy. More... | |
Static Public Member Functions | |
| static vtkPolyhedron * | New () |
| Standard new methods. More... | |
| static vtkTypeBool | IsTypeOf (const char *type) |
| static vtkPolyhedron * | SafeDownCast (vtkObjectBase *o) |
Static Public Member Functions inherited from vtkCell3D | |
| static vtkTypeBool | IsTypeOf (const char *type) |
| static vtkCell3D * | SafeDownCast (vtkObjectBase *o) |
Static Public Member Functions inherited from vtkCell | |
| static vtkTypeBool | IsTypeOf (const char *type) |
| static vtkCell * | SafeDownCast (vtkObjectBase *o) |
Static Public Member Functions inherited from vtkObject | |
| static vtkObject * | New () |
| Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More... | |
| static void | BreakOnError () |
| This method is called when vtkErrorMacro executes. More... | |
| static void | SetGlobalWarningDisplay (int val) |
| This is a global flag that controls whether any debug, warning or error messages are displayed. More... | |
| static void | GlobalWarningDisplayOn () |
| static void | GlobalWarningDisplayOff () |
| static int | GetGlobalWarningDisplay () |
Static Public Member Functions inherited from vtkObjectBase | |
| static vtkTypeBool | IsTypeOf (const char *name) |
| Return 1 if this class type is the same type of (or a subclass of) the named class. More... | |
| static vtkIdType | GetNumberOfGenerationsFromBaseType (const char *name) |
| Given a the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class). More... | |
| static vtkObjectBase * | New () |
| Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More... | |
Protected Member Functions | |
| virtual vtkObjectBase * | NewInstanceInternal () const |
| vtkPolyhedron () | |
| ~vtkPolyhedron () override | |
| int | GenerateEdges () |
| void | GenerateFaces () |
| void | ComputeBounds () |
| void | ComputeParametricCoordinate (const double x[3], double pc[3]) |
| void | ComputePositionFromParametricCoordinate (const double pc[3], double x[3]) |
| void | ConstructPolyData () |
| void | ConstructLocator () |
Protected Member Functions inherited from vtkCell3D | |
| vtkCell3D () | |
| ~vtkCell3D () override | |
Protected Member Functions inherited from vtkCell | |
| vtkCell () | |
| ~vtkCell () override | |
Protected Member Functions inherited from vtkObject | |
| vtkObject () | |
| ~vtkObject () override | |
| void | RegisterInternal (vtkObjectBase *, vtkTypeBool check) override |
| void | UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) override |
| void | InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=nullptr) |
| These methods allow a command to exclusively grab all events. More... | |
| void | InternalReleaseFocus () |
Protected Member Functions inherited from vtkObjectBase | |
| vtkObjectBase () | |
| virtual | ~vtkObjectBase () |
| virtual void | CollectRevisions (ostream &) |
| virtual void | ReportReferences (vtkGarbageCollector *) |
| vtkObjectBase (const vtkObjectBase &) | |
| void | operator= (const vtkObjectBase &) |
Protected Attributes | |
| vtkLine * | Line |
| vtkTriangle * | Triangle |
| vtkQuad * | Quad |
| vtkPolygon * | Polygon |
| vtkTetra * | Tetra |
| vtkIdTypeArray * | GlobalFaces |
| vtkIdTypeArray * | FaceLocations |
| vtkPointIdMap * | PointIdMap |
| int | EdgesGenerated |
| vtkEdgeTable * | EdgeTable |
| vtkIdTypeArray * | Edges |
| vtkIdTypeArray * | EdgeFaces |
| vtkIdTypeArray * | Faces |
| int | FacesGenerated |
| int | BoundsComputed |
| int | PolyDataConstructed |
| vtkPolyData * | PolyData |
| vtkCellArray * | Polys |
| int | LocatorConstructed |
| vtkCellLocator * | CellLocator |
| vtkIdList * | CellIds |
| vtkGenericCell * | Cell |
Protected Attributes inherited from vtkCell3D | |
| vtkOrderedTriangulator * | Triangulator |
| double | MergeTolerance |
| vtkTetra * | ClipTetra |
| vtkDoubleArray * | ClipScalars |
Protected Attributes inherited from vtkCell | |
| double | Bounds [6] |
Protected Attributes inherited from vtkObject | |
| bool | Debug |
| vtkTimeStamp | MTime |
| vtkSubjectHelper * | SubjectHelper |
Protected Attributes inherited from vtkObjectBase | |
| std::atomic< int32_t > | ReferenceCount |
| vtkWeakPointerBase ** | WeakPointers |
Additional Inherited Members | |
Public Attributes inherited from vtkCell | |
| vtkPoints * | Points |
| vtkIdList * | PointIds |
a 3D cell defined by a set of polygonal faces
vtkPolyhedron is a concrete implementation that represents a 3D cell defined by a set of polygonal faces. The polyhedron should be watertight, non-self-intersecting and manifold (each edge is used twice).
Interpolation functions and weights are defined / computed using the method of Mean Value Coordinates (MVC). See the VTK class vtkMeanValueCoordinatesInterpolator for more information.
The class does not require the polyhedron to be convex. However, the polygonal faces must be planar. Non-planar polygonal faces will definitely cause problems, especially in severely warped situations.
Definition at line 57 of file vtkPolyhedron.h.
| typedef vtkCell3D vtkPolyhedron::Superclass |
Definition at line 65 of file vtkPolyhedron.h.
|
protected |
|
overrideprotected |
|
static |
Standard new methods.
|
static |
|
virtual |
Return 1 if this class is the same type of (or a subclass of) the named class.
Returns 0 otherwise. This method works in combination with vtkTypeMacro found in vtkSetGet.h.
Reimplemented from vtkCell3D.
|
static |
|
protectedvirtual |
Reimplemented from vtkCell3D.
| vtkPolyhedron* vtkPolyhedron::NewInstance | ( | ) | const |
|
overridevirtual |
|
inlineoverride |
See vtkCell3D API for description of these methods.
Definition at line 74 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 83 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 84 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 94 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 95 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 100 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 106 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 112 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 118 of file vtkPolyhedron.h.
|
inlineoverride |
Definition at line 124 of file vtkPolyhedron.h.
|
overridevirtual |
|
inlineoverridevirtual |
See the vtkCell API for descriptions of these methods.
Implements vtkCell.
Definition at line 139 of file vtkPolyhedron.h.
|
inlineoverridevirtual |
This cell requires that it be initialized prior to access.
Reimplemented from vtkCell.
Definition at line 144 of file vtkPolyhedron.h.
|
overridevirtual |
Reimplemented from vtkCell.
|
overridevirtual |
A polyhedron is represented internally by a set of polygonal faces.
These faces can be processed to explicitly determine edges.
Implements vtkCell.
|
overridevirtual |
Return the edge cell from the edgeId of the cell.
Implements vtkCell.
|
overridevirtual |
Return the number of faces in the cell.
Implements vtkCell.
|
overridevirtual |
Return the face cell from the faceId of the cell.
Implements vtkCell.
|
overridevirtual |
|
overridevirtual |
Satisfy the vtkCell API.
This method clips the input polyhedron and outputs a new polyhedron. The face information of the output polyhedron is encoded in the output vtkCellArray using a special format: CellLength [nCellFaces, nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...]. Use the static method vtkUnstructuredGrid::DecomposePolyhedronCellArray to convert it into a standard format. Note: the algorithm assumes water-tight polyhedron cells.
Reimplemented from vtkCell3D.
|
overridevirtual |
Satisfy the vtkCell API.
The subId is ignored and zero is always returned. The parametric coordinates pcoords are normalized values in the bounding box of the polyhedron. The weights are determined by evaluating the MVC coordinates. The dist is always zero if the point x[3] is inside the polyhedron; otherwise it's the distance to the surface.
Implements vtkCell.
|
overridevirtual |
|
overridevirtual |
Intersect the line (p1,p2) with a given tolerance tol to determine a point of intersection x[3] with parametric coordinate t along the line.
The parametric coordinates are returned as well (subId can be ignored). Returns true if the line intersects a face.
Implements vtkCell.
Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh.
This method works well for a convex polyhedron but may return wrong result in a concave case. Once triangulation has been performed, the results are saved in ptIds and pts. The ptIds is a vtkIdList with 4xn number of ids (n is the number of result tetrahedrons). The first 4 represent the point ids of the first tetrahedron, the second 4 represents the point ids of the second tetrahedron and so on. The point ids represent global dataset ids. The points of result tetrahedons are stored in pts. Note that there are 4xm output points (m is the number of points in the original polyhedron). A point may be stored multiple times when it is shared by more than one tetrahedrons. The points stored in pts are ordered the same as they are listed in ptIds.
Implements vtkCell.
|
overridevirtual |
Computes derivatives at the point specified by the parameter coordinate.
Current implementation uses all vertices and subId is not used. To accelerate the speed, the future implementation can triangulate and extract the local tetrahedron from subId and pcoords, then evaluate derivatives on the local tetrahedron.
Implements vtkCell.
|
overridevirtual |
Find the boundary face closest to the point defined by the pcoords[3] and subId of the cell (subId can be ignored).
Implements vtkCell.
|
inlineoverridevirtual |
Return the center of the cell in parametric coordinates.
In this cell, the center of the bounding box is returned.
Reimplemented from vtkCell.
Definition at line 353 of file vtkPolyhedron.h.
|
inlineoverridevirtual |
A polyhedron is a full-fledged primary cell.
Reimplemented from vtkCell.
Definition at line 247 of file vtkPolyhedron.h.
|
override |
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
Here we use the MVC calculation process to compute the interpolation functions.
|
override |
|
inlineoverridevirtual |
Methods supporting the definition of faces.
Note that the GetFaces() returns a list of faces in vtkCellArray form; use the method GetNumberOfFaces() to determine the number of faces in the list. The SetFaces() method is also in vtkCellArray form, except that it begins with a leading count indicating the total number of faces in the list.
Reimplemented from vtkCell.
Definition at line 268 of file vtkPolyhedron.h.
|
override |
| int vtkPolyhedron::IsInside | ( | const double | x[3], |
| double | tolerance | ||
| ) |
A method particular to vtkPolyhedron.
It determines whether a point x[3] is inside the polyhedron or not (returns 1 is the point is inside, 0 otherwise). The tolerance is expressed in normalized space; i.e., a fraction of the size of the bounding box.
| bool vtkPolyhedron::IsConvex | ( | ) |
Determine whether or not a polyhedron is convex.
This method is adapted from Devillers et al., "Checking the Convexity of Polytopes and the Planarity of Subdivisions", Computational Geometry, Volume 11, Issues 3 - 4, December 1998, Pages 187 - 208.
| vtkPolyData* vtkPolyhedron::GetPolyData | ( | ) |
Construct polydata if no one exist, then return this->PolyData.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Definition at line 299 of file vtkPolyhedron.h.
|
protected |
Definition at line 300 of file vtkPolyhedron.h.
|
protected |
Definition at line 301 of file vtkPolyhedron.h.
|
protected |
Definition at line 302 of file vtkPolyhedron.h.
|
protected |
Definition at line 303 of file vtkPolyhedron.h.
|
protected |
Definition at line 304 of file vtkPolyhedron.h.
|
protected |
Definition at line 305 of file vtkPolyhedron.h.
|
protected |
Definition at line 312 of file vtkPolyhedron.h.
|
protected |
Definition at line 316 of file vtkPolyhedron.h.
|
protected |
Definition at line 317 of file vtkPolyhedron.h.
|
protected |
Definition at line 318 of file vtkPolyhedron.h.
|
protected |
Definition at line 319 of file vtkPolyhedron.h.
|
protected |
Definition at line 326 of file vtkPolyhedron.h.
|
protected |
Definition at line 327 of file vtkPolyhedron.h.
|
protected |
Definition at line 331 of file vtkPolyhedron.h.
|
protected |
Definition at line 337 of file vtkPolyhedron.h.
|
protected |
Definition at line 338 of file vtkPolyhedron.h.
|
protected |
Definition at line 339 of file vtkPolyhedron.h.
|
protected |
Definition at line 341 of file vtkPolyhedron.h.
|
protected |
Definition at line 342 of file vtkPolyhedron.h.
|
protected |
Definition at line 344 of file vtkPolyhedron.h.
|
protected |
Definition at line 345 of file vtkPolyhedron.h.
1.8.17