VTK
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
dox
Filtering
vtkUnstructuredGrid.h
Go to the documentation of this file.
1
/*=========================================================================
2
3
Program: Visualization Toolkit
4
Module: vtkUnstructuredGrid.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
=========================================================================*/
32
#ifndef __vtkUnstructuredGrid_h
33
#define __vtkUnstructuredGrid_h
34
35
#include "
vtkPointSet.h
"
36
37
class
vtkCellArray
;
38
class
vtkCellLinks
;
39
class
vtkConvexPointSet
;
40
class
vtkEmptyCell
;
41
class
vtkHexahedron
;
42
class
vtkIdList
;
43
class
vtkIdTypeArray
;
44
class
vtkLine
;
45
class
vtkPixel
;
46
class
vtkPolyLine
;
47
class
vtkPolyVertex
;
48
class
vtkPolygon
;
49
class
vtkPyramid
;
50
class
vtkPentagonalPrism
;
51
class
vtkHexagonalPrism
;
52
class
vtkQuad
;
53
class
vtkQuadraticEdge
;
54
class
vtkQuadraticHexahedron
;
55
class
vtkQuadraticWedge
;
56
class
vtkQuadraticPyramid
;
57
class
vtkQuadraticQuad
;
58
class
vtkQuadraticTetra
;
59
class
vtkQuadraticTriangle
;
60
class
vtkTetra
;
61
class
vtkTriangle
;
62
class
vtkTriangleStrip
;
63
class
vtkUnsignedCharArray
;
64
class
vtkVertex
;
65
class
vtkVoxel
;
66
class
vtkWedge
;
67
class
vtkTriQuadraticHexahedron
;
68
class
vtkQuadraticLinearWedge
;
69
class
vtkQuadraticLinearQuad
;
70
class
vtkBiQuadraticQuad
;
71
class
vtkBiQuadraticQuadraticWedge
;
72
class
vtkBiQuadraticQuadraticHexahedron
;
73
class
vtkBiQuadraticTriangle
;
74
class
vtkCubicLine
;
75
class
vtkPolyhedron
;
76
class
vtkIdTypeArray
;
77
78
class
VTK_FILTERING_EXPORT
vtkUnstructuredGrid
:
public
vtkPointSet
79
{
80
public
:
81
static
vtkUnstructuredGrid
*
New
();
82
83
vtkTypeMacro(
vtkUnstructuredGrid
,
vtkPointSet
);
84
void
PrintSelf
(ostream& os,
vtkIndent
indent);
85
87
88
int
GetDataObjectType
() {
return
VTK_UNSTRUCTURED_GRID
;};
89
virtual
void
Allocate(
vtkIdType
numCells=1000,
int
extSize=1000);
91
100
vtkIdType
InsertNextCell(
int
type
,
vtkIdType
npts,
vtkIdType
*ptIds);
101
108
vtkIdType
InsertNextCell(
int
type
,
vtkIdList
*ptIds);
109
110
// Desciption:
111
// Insert/create a polyhedron cell. npts is the number of unique points in
112
// the cell. pts is the list of the unique cell point Ids. nfaces is the
113
// number of faces in the cell. faces is the face-stream
114
// [numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...].
115
// All point Ids are global.
116
vtkIdType
InsertNextCell(
int
type
,
vtkIdType
npts,
vtkIdType
*ptIds,
117
vtkIdType
nfaces,
vtkIdType
*faces);
118
120
121
void
Reset();
122
virtual
void
CopyStructure
(
vtkDataSet
*ds);
123
vtkIdType
GetNumberOfCells
();
124
virtual
vtkCell
*
GetCell
(
vtkIdType
cellId);
125
virtual
void
GetCell
(
vtkIdType
cellId,
vtkGenericCell
*cell);
126
virtual
void
GetCellBounds
(
vtkIdType
cellId,
double
bounds[6]);
127
virtual
void
GetCellPoints
(
vtkIdType
cellId,
vtkIdList
*ptIds);
128
void
GetPointCells
(
vtkIdType
ptId,
vtkIdList
*cellIds);
130
131
int
GetCellType
(
vtkIdType
cellId);
132
vtkUnsignedCharArray
*
GetCellTypesArray
() {
return
this->Types; }
133
vtkIdTypeArray
*
GetCellLocationsArray
() {
return
this->Locations; }
134
void
Squeeze
();
135
void
Initialize
();
136
int
GetMaxCellSize
();
137
void
BuildLinks();
138
vtkCellLinks
*
GetCellLinks
() {
return
this->Links;};
139
virtual
void
GetCellPoints
(
vtkIdType
cellId,
vtkIdType
& npts,
140
vtkIdType
* &pts);
141
147
void
GetFaceStream(
vtkIdType
cellId,
vtkIdList
*ptIds);
148
154
void
GetFaceStream(
vtkIdType
cellId,
vtkIdType
& nfaces,
vtkIdType
* &ptIds);
155
157
169
void
SetCells(
int
type
,
vtkCellArray
*cells);
170
void
SetCells(
int
*types,
vtkCellArray
*cells);
171
void
SetCells(
vtkUnsignedCharArray
*cellTypes,
vtkIdTypeArray
*cellLocations,
172
vtkCellArray
*cells);
173
void
SetCells(
vtkUnsignedCharArray
*cellTypes,
vtkIdTypeArray
*cellLocations,
174
vtkCellArray
*cells,
vtkIdTypeArray
*faceLocations,
175
vtkIdTypeArray
*faces);
177
178
vtkCellArray
*
GetCells
() {
return
this->Connectivity;};
179
void
ReplaceCell(
vtkIdType
cellId,
int
npts,
vtkIdType
*pts);
180
vtkIdType
InsertNextLinkedCell(
int
type
,
int
npts,
vtkIdType
*pts);
181
void
RemoveReferenceToCell(
vtkIdType
ptId,
vtkIdType
cellId);
182
void
AddReferenceToCell(
vtkIdType
ptId,
vtkIdType
cellId);
183
void
ResizeCellList(
vtkIdType
ptId,
int
size
);
184
186
189
virtual
void
GetCellNeighbors
(
vtkIdType
cellId,
vtkIdList
*ptIds,
190
vtkIdList
*cellIds);
192
196
void
GetUpdateExtent
(
int
&piece,
int
&numPieces,
int
&ghostLevel);
197
199
200
virtual
int
*
GetUpdateExtent
();
201
virtual
void
GetUpdateExtent
(
int
& x0,
int
& x1,
int
& y0,
int
& y1,
202
int
& z0,
int
& z1);
203
virtual
void
GetUpdateExtent
(
int
extent
[6]);
205
207
209
virtual
int
GetPiece();
210
virtual
int
GetNumberOfPieces();
212
214
virtual
int
GetGhostLevel();
215
221
unsigned
long
GetActualMemorySize
();
222
224
225
virtual
void
ShallowCopy
(
vtkDataObject
*src);
226
virtual
void
DeepCopy
(
vtkDataObject
*src);
228
232
void
GetIdsOfCellsOfType(
int
type
,
vtkIdTypeArray
*array);
233
235
int
IsHomogeneous();
236
239
void
RemoveGhostCells(
int
level
);
240
241
//BTX
243
244
static
vtkUnstructuredGrid
*
GetData
(
vtkInformation
*
info
);
245
static
vtkUnstructuredGrid
*
GetData
(
vtkInformationVector
* v,
int
i=0);
246
//ETX
248
250
vtkIdType
*GetFaces(
vtkIdType
cellId);
251
253
254
vtkIdTypeArray
*
GetFaces
(){
return
this->Faces;};
255
vtkIdTypeArray
*
GetFaceLocations
(){
return
this->FaceLocations;};
257
264
int
InitializeFacesRepresentation(
vtkIdType
numPrevCells);
265
267
276
static
void
DecomposeAPolyhedronCell(
vtkCellArray
*polyhedronCellArray,
277
vtkIdType
& nCellpts,
278
vtkIdType
& nCellfaces,
279
vtkCellArray
*cellArray,
280
vtkIdTypeArray
*faces);
282
283
static
void
DecomposeAPolyhedronCell(
vtkIdType
* polyhedronCellStream,
284
vtkIdType
& nCellpts,
285
vtkIdType
& nCellfaces,
286
vtkCellArray
*cellArray,
287
vtkIdTypeArray
*faces);
288
290
299
static
void
DecomposeAPolyhedronCell(
vtkIdType
nCellFaces,
300
vtkIdType
* inFaceStream,
301
vtkIdType
& nCellpts,
302
vtkCellArray
* cellArray,
303
vtkIdTypeArray
* faces);
305
307
311
static
void
ConvertFaceStreamPointIds(
vtkIdList
* faceStream,
312
vtkIdType
* idMap);
314
316
320
static
void
ConvertFaceStreamPointIds(
vtkIdType
nfaces,
321
vtkIdType
* faceStream,
322
vtkIdType
* idMap);
324
325
326
protected
:
327
vtkUnstructuredGrid
();
328
~
vtkUnstructuredGrid
();
329
330
// used by GetCell method
331
vtkVertex
*
Vertex
;
332
vtkPolyVertex
*
PolyVertex
;
333
vtkLine
*
Line
;
334
vtkPolyLine
*
PolyLine
;
335
vtkTriangle
*
Triangle
;
336
vtkTriangleStrip
*
TriangleStrip
;
337
vtkPixel
*
Pixel
;
338
vtkQuad
*
Quad
;
339
vtkPolygon
*
Polygon
;
340
vtkTetra
*
Tetra
;
341
vtkVoxel
*
Voxel
;
342
vtkHexahedron
*
Hexahedron
;
343
vtkWedge
*
Wedge
;
344
vtkPyramid
*
Pyramid
;
345
vtkPentagonalPrism
*
PentagonalPrism
;
346
vtkHexagonalPrism
*
HexagonalPrism
;
347
vtkQuadraticEdge
*
QuadraticEdge
;
348
vtkQuadraticTriangle
*
QuadraticTriangle
;
349
vtkQuadraticQuad
*
QuadraticQuad
;
350
vtkQuadraticTetra
*
QuadraticTetra
;
351
vtkQuadraticHexahedron
*
QuadraticHexahedron
;
352
vtkQuadraticWedge
*
QuadraticWedge
;
353
vtkQuadraticPyramid
*
QuadraticPyramid
;
354
vtkQuadraticLinearQuad
*
QuadraticLinearQuad
;
355
vtkBiQuadraticQuad
*
BiQuadraticQuad
;
356
vtkTriQuadraticHexahedron
*
TriQuadraticHexahedron
;
357
vtkQuadraticLinearWedge
*
QuadraticLinearWedge
;
358
vtkBiQuadraticQuadraticWedge
*
BiQuadraticQuadraticWedge
;
359
vtkBiQuadraticQuadraticHexahedron
*
BiQuadraticQuadraticHexahedron
;
360
vtkBiQuadraticTriangle
*
BiQuadraticTriangle
;
361
vtkCubicLine
*
CubicLine
;
362
vtkConvexPointSet
*
ConvexPointSet
;
363
vtkPolyhedron
*
Polyhedron
;
364
vtkEmptyCell
*
EmptyCell
;
365
366
// points inherited
367
// point data (i.e., scalars, vectors, normals, tcoords) inherited
368
vtkCellArray
*
Connectivity
;
369
vtkCellLinks
*
Links
;
370
vtkUnsignedCharArray
*
Types
;
371
vtkIdTypeArray
*
Locations
;
372
373
// Special support for polyhedra/cells with explicit face representations.
374
// The Faces class represents polygonal faces using a modified vtkCellArray
375
// structure. Each cell face list begins with the total number of faces in
376
// the cell, followed by a vtkCellArray data organization
377
// (n,i,j,k,n,i,j,k,...).
378
vtkIdTypeArray
*
Faces
;
379
vtkIdTypeArray
*
FaceLocations
;
380
381
private
:
382
// Hide these from the user and the compiler.
383
vtkUnstructuredGrid
(
const
vtkUnstructuredGrid
&);
// Not implemented.
384
void
operator=(
const
vtkUnstructuredGrid
&);
// Not implemented.
385
386
void
Cleanup();
387
389
390
VTK_LEGACY(
void
GetCellNeighbors
(
vtkIdType
cellId,
vtkIdList
& ptIds,
vtkIdList
& cellIds));
391
};
393
394
#endif
Generated on Sun Sep 9 2012 13:03:29 for VTK by
1.8.1.2