VTK
vtkMappedUnstructuredGrid.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkMappedUnstructuredGrid.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 =========================================================================*/
139 #ifndef vtkMappedUnstructuredGrid_h
140 #define vtkMappedUnstructuredGrid_h
141 
142 #include "vtkUnstructuredGridBase.h"
143 
144 #include "vtkMappedUnstructuredGridCellIterator.h" // For default cell iterator
145 #include "vtkNew.h" // For vtkNew
146 #include "vtkSmartPointer.h" // For vtkSmartPointer
147 
148 template <class Implementation,
152 {
154 public:
156  typedef Implementation ImplementationType;
157  typedef CellIterator CellIteratorType;
158 
159  // Virtuals from various base classes:
160  void PrintSelf(ostream &os, vtkIndent indent);
161  void CopyStructure(vtkDataSet *pd);
162  void ShallowCopy(vtkDataObject *src);
164  vtkCell* GetCell(vtkIdType cellId);
165  void GetCell(vtkIdType cellId, vtkGenericCell *cell);
166  int GetCellType(vtkIdType cellId);
167  void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds);
169  void GetPointCells(vtkIdType ptId, vtkIdList *cellIds);
170  int GetMaxCellSize();
171  void GetIdsOfCellsOfType(int type, vtkIdTypeArray *array);
172  int IsHomogeneous();
173  void Allocate(vtkIdType numCells, int extSize = 1000);
174  vtkIdType InsertNextCell(int type, vtkIdList *ptIds);
175  vtkIdType InsertNextCell(int type, vtkIdType npts, vtkIdType *ptIds);
176  vtkIdType InsertNextCell(int type, vtkIdType npts, vtkIdType *ptIds,
177  vtkIdType nfaces, vtkIdType *faces);
178  void ReplaceCell(vtkIdType cellId, int npts, vtkIdType *pts);
179  vtkMTimeType GetMTime() VTK_OVERRIDE;
180 
181  void SetImplementation(ImplementationType *impl);
182  ImplementationType *GetImplementation();
183 
184 protected:
186  ~vtkMappedUnstructuredGrid();
187 
188  // For convenience...
189  typedef vtkMappedUnstructuredGrid<Implementation, CellIterator> ThisType;
190 
191  vtkSmartPointer<ImplementationType> Impl;
192 
193 private:
194  vtkMappedUnstructuredGrid(const vtkMappedUnstructuredGrid &) VTK_DELETE_FUNCTION;
195  void operator=(const vtkMappedUnstructuredGrid &) VTK_DELETE_FUNCTION;
196 
197  vtkNew<vtkGenericCell> TempCell;
198 };
199 
200 #include "vtkMappedUnstructuredGrid.txx"
201 
202 // We need to fake the superclass for the wrappers, otherwise they will choke on
203 // the template:
204 #ifndef __VTK_WRAP__
205 
206 #define vtkMakeExportedMappedUnstructuredGrid(_className, _impl, _exportDecl) \
207 class _exportDecl _className : \
208  public vtkMappedUnstructuredGrid<_impl> \
209 { \
210 public: \
211  vtkTypeMacro(_className, \
212  vtkMappedUnstructuredGrid<_impl>) \
213  static _className* New(); \
214 protected: \
215  _className() \
216  { \
217  _impl *i = _impl::New(); \
218  this->SetImplementation(i); \
219  i->Delete(); \
220  } \
221  ~_className() {} \
222 private: \
223  _className(const _className&); \
224  void operator=(const _className&); \
225 };
226 
227 #define vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl) \
228 class _exportDecl _className : \
229  public vtkMappedUnstructuredGrid<_impl, _cIter> \
230 { \
231 public: \
232  vtkTypeMacro(_className, \
233  vtkMappedUnstructuredGrid<_impl, _cIter>) \
234  static _className* New(); \
235 protected: \
236  _className() \
237  { \
238  _impl *i = _impl::New(); \
239  this->SetImplementation(i); \
240  i->Delete(); \
241  } \
242  ~_className() {} \
243 private: \
244  _className(const _className&); \
245  void operator=(const _className&); \
246 };
247 
248 #else // __VTK_WRAP__
249 
250 #define vtkMakeExportedMappedUnstructuredGrid(_className, _impl, _exportDecl) \
251  class _exportDecl _className : \
252  public vtkUnstructuredGridBase \
253  { \
254 public: \
255  vtkTypeMacro(_className, vtkUnstructuredGridBase) \
256  static _className* New(); \
257 protected: \
258  _className() {} \
259  ~_className() {} \
260 private: \
261  _className(const _className&); \
262  void operator=(const _className&); \
263  };
264 
265 #define vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl) \
266  class _exportDecl _className : \
267  public vtkUnstructuredGridBase \
268  { \
269 public: \
270  vtkTypeMacro(_className, vtkUnstructuredGridBase) \
271  static _className* New(); \
272 protected: \
273  _className() {} \
274  ~_className() {} \
275 private: \
276  _className(const _className&); \
277  void operator=(const _className&); \
278  };
279 
280 #endif // __VTK_WRAP__
281 
282 #define vtkMakeMappedUnstructuredGrid(_className, _impl) \
283  vtkMakeExportedMappedUnstructuredGrid(_className, _impl, )
284 
285 #define vtkMakeMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl) \
286  vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, )
287 
288 #endif //vtkMappedUnstructuredGrid_h
289 
290 // VTK-HeaderTest-Exclude: vtkMappedUnstructuredGrid.h
void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)
Topological inquiry to get points defining cell.
void CopyStructure(vtkDataSet *pd)
Copy the geometric and topological structure of an object.
Default cell iterator for vtkMappedUnstructuredGrid.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
vtkTemplateTypeMacro(SelfType, vtkUnstructuredGridBase) typedef Implementation ImplementationType
void GetIdsOfCellsOfType(int type, vtkIdTypeArray *array)
Fill vtkIdTypeArray container with list of cell Ids.
dynamic, self-adjusting array of vtkIdType
Hold a reference to a vtkObjectBase instance.
int GetMaxCellSize()
Convenience method returns largest cell size in dataset.
int vtkIdType
Definition: vtkType.h:287
void SetImplementation(ImplementationType *impl)
vtkCellIterator * NewCellIterator()
Return an iterator that traverses the cells in this data set.
vtkIdType InsertNextCell(int type, vtkIdList *ptIds)
Insert/create cell in object by a list of point ids defining cell topology.
vtkTypeUInt64 vtkMTimeType
Definition: vtkType.h:248
provides thread-safe access to cells
vtkMTimeType GetMTime() override
Datasets are composite objects and need to check each part for MTime THIS METHOD IS THREAD SAFE...
void ReplaceCell(vtkIdType cellId, int npts, vtkIdType *pts)
Replace the points defining cell "cellId" with a new set of points.
abstract class to specify cell behavior
Definition: vtkCell.h:59
a simple class to control print indentation
Definition: vtkIndent.h:39
ImplementationType * GetImplementation()
list of point or cell ids
Definition: vtkIdList.h:36
vtkIdType GetNumberOfCells()
Determine the number of cells composing the dataset.
int GetCellType(vtkIdType cellId)
Get type of cell with cellId such that: 0 <= cellId < NumberOfCells.
void GetPointCells(vtkIdType ptId, vtkIdList *cellIds)
Topological inquiry to get cells using point.
Allocate and hold a VTK object.
Definition: vtkNew.h:67
int IsHomogeneous()
Traverse cells and determine if cells are all of the same type.
void ShallowCopy(vtkDataObject *src)
Shallow and Deep copy.
vtkSmartPointer< ImplementationType > Impl
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
Efficient cell iterator for vtkDataSet topologies.
vtkCell * GetCell(vtkIdType cellId)
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
general representation of visualization data
Definition: vtkDataObject.h:64
dataset represents arbitrary combinations of all possible cell types.
void Allocate(vtkIdType numCells, int extSize=1000)
Allocate memory for the number of cells indicated.
Allows datasets with arbitrary storage layouts to be used with VTK.