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 =========================================================================*/
15 /**
16  * @class   vtkMappedUnstructuredGrid
17  * @brief   Allows datasets with arbitrary storage
18  * layouts to be used with VTK.
19  *
20  *
21  * This class fulfills the vtkUnstructuredGridBase API while delegating to a
22  * arbitrary implementation of the dataset topology. The purpose of
23  * vtkMappedUnstructuredGrid is to allow external data structures to be used
24  * directly in a VTK pipeline, e.g. for in-situ analysis of a running
25  * simulation.
26  *
27  * When introducing an external data structure into VTK, there are 3 principle
28  * components of the dataset to consider:
29  * - Points
30  * - Cells (topology)
31  * - Point/Cell attributes
32  *
33  * Points and attributes can be handled by subclassing vtkMappedDataArray and
34  * implementing that interface to adapt the external data structures into VTK.
35  * The vtkMappedDataArray subclasses can then be used as the vtkPoints's Data
36  * member (for points/nodes) or added directly to vtkPointData, vtkCellData, or
37  * vtkFieldData for attribute information. Filters used in the pipeline will
38  * need to be modified to remove calls to vtkDataArray::GetVoidPointer and use
39  * a suitable vtkArrayDispatch instead.
40  *
41  * Introducing an arbitrary topology implementation into VTK requires the use of
42  * the vtkMappedUnstructuredGrid class. Unlike the data array counterpart, the
43  * mapped unstructured grid is not subclassed, but rather takes an adaptor
44  * class as a template argument. This is to allow cheap shallow copies of the
45  * data by passing the reference-counted implementation object to new instances
46  * of vtkMappedUnstructuredGrid.
47  *
48  * The implementation class should derive from vtkObject (for reference
49  * counting) and implement the usual vtkObject API requirements, such as a
50  * static New() method and PrintSelf function. The following methods must also
51  * be implemented:
52  * - vtkIdType GetNumberOfCells()
53  * - int GetCellType(vtkIdType cellId)
54  * - void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds)
55  * - void GetPointCells(vtkIdType ptId, vtkIdList *cellIds)
56  * - int GetMaxCellSize()
57  * - void GetIdsOfCellsOfType(int type, vtkIdTypeArray *array)
58  * - int IsHomogeneous()
59  * - void Allocate(vtkIdType numCells, int extSize = 1000)
60  * - vtkIdType InsertNextCell(int type, vtkIdList *ptIds)
61  * - vtkIdType InsertNextCell(int type, vtkIdType npts, const vtkIdType ptIds[])
62  * - vtkIdType InsertNextCell(int type, vtkIdType npts, const vtkIdType ptIds[],
63  *                            vtkIdType nfaces, const vtkIdType faces[])
64  * - void ReplaceCell(vtkIdType cellId, int npts, const vtkIdType pts[])
65  *
66  * These methods should provide the same functionality as defined in
67  * vtkUnstructuredGrid. See that class's documentation for more information.
68  *
69  * Note that since the implementation class is used as a compile-time template
70  * parameter in vtkMappedUnstructuredGrid, the above methods do not need be
71  * virtuals. The compiler will statically bind the calls, making dynamic vtable
72  * lookups unnecessary and giving a slight performance boost.
73  *
74  * Adapting a filter or algorithm to safely traverse the
75  * vtkMappedUnstructuredGrid's topology requires removing calls the following
76  * implementation-dependent vtkUnstructuredGrid methods:
77  * - vtkUnstructuredGrid::GetCellTypesArray()
78  * - vtkUnstructuredGrid::GetCellLocationsArray()
79  * - vtkUnstructuredGrid::GetCellLinks()
80  * - vtkUnstructuredGrid::GetCells()
81  * Access to the values returned by these methods should be replaced by the
82  * equivalent random-access lookup methods in the vtkUnstructuredGridBase API,
83  * or use vtkCellIterator (see vtkDataSet::NewCellIterator) for sequential
84  * access.
85  *
86  * A custom vtkCellIterator implementation may be specified for a particular
87  * vtkMappedUnstructuredGrid as the second template parameter. By default,
88  * vtkMappedUnstructuredGridCellIterator will be used, which increments an
89  * internal cell id counter and performs random-access lookup as-needed. More
90  * efficient implementations may be used with data structures better suited for
91  * sequential access, see vtkUnstructuredGridCellIterator for an example.
92  *
93  * A set of four macros are provided to generate a concrete subclass of
94  * vtkMappedUnstructuredGrid with a specified implementation, cell iterator,
95  * and export declaration. They are:
96  * - vtkMakeMappedUnstructuredGrid(_className, _impl)
97  *   - Create a subclass of vtkMappedUnstructuredGrid using _impl implementation
98  *     that is named _className.
99  * - vtkMakeMappedUnstructuredGridWithIter(_className, _impl, _cIter)
100  *   - Create a subclass of vtkMappedUnstructuredGrid using _impl implementation
101  *     and _cIter vtkCellIterator that is named _className.
102  * - vtkMakeExportedMappedUnstructuredGrid(_className, _impl, _exportDecl)
103  *   - Create a subclass of vtkMappedUnstructuredGrid using _impl implementation
104  *     that is named _className. _exportDecl is used to decorate the class
105  *     declaration.
106  * - vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl)
107  *   - Create a subclass of vtkMappedUnstructuredGrid using _impl implementation
108  *     and _cIter vtkCellIterator that is named _className. _exportDecl is used
109  *     to decorate the class declaration.
110  *
111  * To instantiate a vtkMappedUnstructuredGrid subclass created by the above
112  * macro, the follow pattern is encouraged:
113  *
114  * @code
115  * MyGrid.h:
116  * ----------------------------------------------------------------------
117  * class MyGridImplementation : public vtkObject
118  * {
119  * public:
120  *   ... (vtkObject required API) ...
121  *   ... (vtkMappedUnstructuredGrid Implementation required API) ...
122  *   void SetImplementationDetails(...raw data from external source...);
123  * };
124  *
125  * vtkMakeMappedUnstructuredGrid(MyGrid, MyGridImplementation)
126  *
127  * SomeSource.cxx
128  * ----------------------------------------------------------------------
129  * vtkNew<MyGrid> grid;
130  * grid->GetImplementation()->SetImplementationDetails(...);
131  * // grid is now ready to use.
132  * @endcode
133  *
134  * The vtkCPExodusIIElementBlock class provides an example of
135  * vtkMappedUnstructuredGrid usage, adapting the Exodus II data structures for
136  * the VTK pipeline.
137 */
138 
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,
149           class CellIterator = vtkMappedUnstructuredGridCellIterator<Implementation> >
150 class vtkMappedUnstructuredGrid:
151     public vtkUnstructuredGridBase
152 {
153   typedef vtkMappedUnstructuredGrid<Implementation, CellIterator> SelfType;
154 public:
155   vtkTemplateTypeMacro(SelfType, vtkUnstructuredGridBase)
156   typedef Implementation ImplementationType;
157   typedef CellIterator CellIteratorType;
158 
159   // Virtuals from various base classes:
160   void PrintSelf(ostream &os, vtkIndent indent) override;
161   void CopyStructure(vtkDataSet *pd) override;
162   void ShallowCopy(vtkDataObject *src) override;
163   vtkIdType GetNumberOfCells() override;
164   using vtkDataSet::GetCell;
165   vtkCell* GetCell(vtkIdType cellId) override;
166   void GetCell(vtkIdType cellId, vtkGenericCell *cell) override;
167   int GetCellType(vtkIdType cellId) override;
168   void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds) override;
169   vtkCellIterator* NewCellIterator() override;
170   void GetPointCells(vtkIdType ptId, vtkIdList *cellIds) override;
171   int GetMaxCellSize() override;
172   void GetIdsOfCellsOfType(int type, vtkIdTypeArray *array) override;
173   int IsHomogeneous() override;
174   void Allocate(vtkIdType numCells, int extSize = 1000) override;
175   vtkMTimeType GetMTime() override;
176 
177   void SetImplementation(ImplementationType *impl);
178   ImplementationType *GetImplementation();
179 
180 protected:
181   vtkMappedUnstructuredGrid();
182   ~vtkMappedUnstructuredGrid() override;
183 
184   // For convenience...
185   typedef vtkMappedUnstructuredGrid<Implementation, CellIterator> ThisType;
186 
187   vtkSmartPointer<ImplementationType> Impl;
188 
189   vtkIdType InternalInsertNextCell(int type, vtkIdType npts, const vtkIdType ptIds[]) override;
190   vtkIdType InternalInsertNextCell(int type, vtkIdList *ptIds) override;
191   vtkIdType InternalInsertNextCell(int type, vtkIdType npts, const vtkIdType ptIds[],
192     vtkIdType nfaces, const vtkIdType faces[]) override;
193   void InternalReplaceCell(vtkIdType cellId, int npts, const vtkIdType pts[]) override;
194 
195 private:
196   vtkMappedUnstructuredGrid(const vtkMappedUnstructuredGrid &) = delete;
197   void operator=(const vtkMappedUnstructuredGrid &) = delete;
198 
199   vtkNew<vtkGenericCell> TempCell;
200 };
201 
202 #include "vtkMappedUnstructuredGrid.txx"
203 
204 // We need to fake the superclass for the wrappers, otherwise they will choke on
205 // the template:
206 #ifndef __VTK_WRAP__
207 
208 #define vtkMakeExportedMappedUnstructuredGrid(_className, _impl, _exportDecl) \
209 class _exportDecl _className : \
210     public vtkMappedUnstructuredGrid<_impl> \
211 { \
212 public: \
213   vtkTypeMacro(_className, \
214                vtkMappedUnstructuredGrid<_impl>) \
215   static _className* New(); \
216 protected: \
217   _className() \
218   { \
219     _impl *i = _impl::New(); \
220     this->SetImplementation(i); \
221     i->Delete(); \
222   } \
223   ~_className() override {} \
224 private: \
225   _className(const _className&); \
226   void operator=(const _className&); \
227 };
228 
229 #define vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl) \
230 class _exportDecl _className : \
231   public vtkMappedUnstructuredGrid<_impl, _cIter> \
232 { \
233 public: \
234   vtkTypeMacro(_className, \
235                vtkMappedUnstructuredGrid<_impl, _cIter>) \
236   static _className* New(); \
237 protected: \
238   _className() \
239   { \
240     _impl *i = _impl::New(); \
241     this->SetImplementation(i); \
242     i->Delete(); \
243   } \
244   ~_className() override {} \
245 private: \
246   _className(const _className&); \
247   void operator=(const _className&); \
248 };
249 
250 #else // __VTK_WRAP__
251 
252 #define vtkMakeExportedMappedUnstructuredGrid(_className, _impl, _exportDecl) \
253   class _exportDecl _className : \
254   public vtkUnstructuredGridBase \
255   { \
256 public: \
257   vtkTypeMacro(_className, vtkUnstructuredGridBase) \
258   static _className* New(); \
259 protected: \
260   _className() {} \
261   ~_className() override {} \
262 private: \
263   _className(const _className&); \
264   void operator=(const _className&); \
265   };
266 
267 #define vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl) \
268   class _exportDecl _className : \
269   public vtkUnstructuredGridBase \
270   { \
271 public: \
272   vtkTypeMacro(_className, vtkUnstructuredGridBase) \
273   static _className* New(); \
274 protected: \
275   _className() {} \
276   ~_className() override {} \
277 private: \
278   _className(const _className&); \
279   void operator=(const _className&); \
280   };
281 
282 #endif // __VTK_WRAP__
283 
284 #define vtkMakeMappedUnstructuredGrid(_className, _impl) \
285   vtkMakeExportedMappedUnstructuredGrid(_className, _impl, )
286 
287 #define vtkMakeMappedUnstructuredGridWithIter(_className, _impl, _cIter, _exportDecl) \
288   vtkMakeExportedMappedUnstructuredGridWithIter(_className, _impl, _cIter, )
289 
290 #endif //vtkMappedUnstructuredGrid_h
291 
292 // VTK-HeaderTest-Exclude: vtkMappedUnstructuredGrid.h
293