SENSEI
A frame work for generic in situ analytics
SVTKUtils.h
Go to the documentation of this file.
1 #ifndef SVTKUtils_h
2 #define SVTKUtils_h
3 
4 /// @file
5 
6 #include "MeshMetadata.h"
7 #include "Error.h"
8 
9 class svtkDataSet;
10 class svtkDataObject;
11 class svtkFieldData;
12 class svtkCellData;
13 class svtkPointData;
14 class svtkDataArray;
15 class svtkTypeInt64Array;
16 class svtkTypeInt32Array;
17 class svtkCellArray;
18 class svtkImageData;
19 class svtkUniformGrid;
20 class svtkRectilinearGrid;
21 class svtkStructuredGrid;
22 class svtkPolyData;
23 class svtkUnstructuredGrid;
24 class svtkMultiBlockDataSet;
25 class svtkOverlappingAMR;
26 class svtkDataSetAttributes;
27 class svtkCompositeDataSet;
28 
29 class vtkDataArray;
30 class vtkTypeInt64Array;
31 class vtkTypeInt32Array;
32 class vtkCellArray;
33 class vtkFieldData;
34 class vtkPointData;
35 class vtkCellData;
36 class vtkPoints;
37 class vtkImageData;
38 class vtkUniformGrid;
39 class vtkRectilinearGrid;
40 class vtkStructuredGrid;
41 class vtkPolyData;
42 class vtkUnstructuredGrid;
43 class vtkMultiBlockDataSet;
44 class vtkOverlappingAMR;
45 class vtkDataObject;
46 class vtkDataSet;
47 
48 #include <svtkDataArray.h>
49 #include <svtkAOSDataArrayTemplate.h>
50 #include <svtkSOADataArrayTemplate.h>
51 
52 #include <svtkSmartPointer.h>
53 #include <svtkAOSDataArrayTemplate.h>
54 #include <svtkCellArray.h>
55 
56 #include <functional>
57 #include <vector>
58 #include <mpi.h>
59 
60 #define svtkCellTemplateMacro(code) \
61  svtkTemplateMacroCase(SVTK_LONG_LONG, long long, code); \
62  svtkTemplateMacroCase(SVTK_UNSIGNED_LONG_LONG, unsigned long long, code); \
63  svtkTemplateMacroCase(SVTK_ID_TYPE, svtkIdType, code); \
64  svtkTemplateMacroCase(SVTK_LONG, long, code); \
65  svtkTemplateMacroCase(SVTK_UNSIGNED_LONG, unsigned long, code); \
66  svtkTemplateMacroCase(SVTK_INT, int, code); \
67  svtkTemplateMacroCase(SVTK_UNSIGNED_INT, unsigned int, code);
68 
69 #define svtkTemplateMacroFP(call) \
70  svtkTemplateMacroCase(SVTK_DOUBLE, double, call); \
71  svtkTemplateMacroCase(SVTK_FLOAT, float, call);
72 
73 
74 #if defined(SENSEI_ENABLE_VTK_CORE)
75 
76 #include <vtkSetGet.h>
77 
78 #define vtkCellTemplateMacro(code) \
79  vtkTemplateMacroCase(VTK_LONG_LONG, long long, code); \
80  vtkTemplateMacroCase(VTK_UNSIGNED_LONG_LONG, unsigned long long, code); \
81  vtkTemplateMacroCase(VTK_ID_TYPE, vtkIdType, code); \
82  vtkTemplateMacroCase(VTK_LONG, long, code); \
83  vtkTemplateMacroCase(VTK_UNSIGNED_LONG, unsigned long, code); \
84  vtkTemplateMacroCase(VTK_INT, int, code); \
85  vtkTemplateMacroCase(VTK_UNSIGNED_INT, unsigned int, code);
86 
87 #define vtkTemplateMacroFP(call) \
88  vtkTemplateMacroCase(VTK_DOUBLE, double, call); \
89  vtkTemplateMacroCase(VTK_FLOAT, float, call);
90 #endif
91 
92 
93 using svtkCompositeDataSetPtr = svtkSmartPointer<svtkCompositeDataSet>;
94 
95 namespace sensei
96 {
97 
98 /** A collection of generally useful funcitons implementing common access
99  * patterns or operations on SVTK data structures.
100  */
101 namespace SVTKUtils
102 {
103 /** given a svtkDataArray get a pointer to underlying data
104  * this handles access from SVTK's AOS and SOA layouts. For
105  * SOA layout only single component arrays should be passed.
106  */
107 template <typename SVTK_TT>
109 {
110  using AOS_ARRAY_TT = svtkAOSDataArrayTemplate<SVTK_TT>;
111  using SOA_ARRAY_TT = svtkSOADataArrayTemplate<SVTK_TT>;
112 
113  AOS_ARRAY_TT *aosDa = nullptr;
114  SOA_ARRAY_TT *soaDa = nullptr;
115 
116  if ((aosDa = dynamic_cast<AOS_ARRAY_TT*>(da)))
117  {
118  return aosDa->GetPointer(0);
119  }
120  else if ((soaDa = dynamic_cast<SOA_ARRAY_TT*>(da)))
121  {
122  return soaDa->GetPointer(0);
123  }
124 
125  SENSEI_ERROR("Invalid svtkDataArray "
126  << (da ? da->GetClassName() : "nullptr"))
127  return nullptr;
128 }
129 
130 /// given a SVTK type enum returns the sizeof that type
131 SENSEI_EXPORT
132 unsigned int Size(int svtkt);
133 
134 /// given a SVTK data object enum returns true if it a legacy object
135 SENSEI_EXPORT
136 int IsLegacyDataObject(int code);
137 
138 /// givne a SVTK data object enum constructs an instance
139 SENSEI_EXPORT
140 svtkDataObject *NewDataObject(int code);
141 
142 /** returns the enum value given an association name. where name can be one of:
143  * point, cell or, field
144  */
145 SENSEI_EXPORT
146 int GetAssociation(std::string assocStr, int &assoc);
147 
148 /// returns the name of the association, point, cell or field
149 SENSEI_EXPORT
150 const char *GetAttributesName(int association);
151 
152 /** returns the container for the associations: svtkPointData, svtkCellData, or
153  * svtkFieldData
154  */
155 SENSEI_EXPORT
156 svtkFieldData *GetAttributes(svtkDataSet *dobj, int association);
157 
158 /** A callback that processes input and output datasets.
159  * return 0 for success, > zero to stop without error, < zero to stop with error
160  */
161 using BinaryDatasetFunction = std::function<int(svtkDataSet*, svtkDataSet*)>;
162 
163 /** Applies the function to leaves of the structurally equivalent
164  * input and output data objects.
165  */
166 SENSEI_EXPORT
167 int Apply(svtkDataObject *input, svtkDataObject *output,
168  BinaryDatasetFunction &func);
169 
170 /** A callback that processes input and output datasets
171  * return 0 for success, > zero to stop without error, < zero to stop with error
172  */
173 using DatasetFunction = std::function<int(svtkDataSet*)>;
174 
175 /** Applies the function to the data object. The function is called once for
176  * each leaf dataset
177  */
178 SENSEI_EXPORT
179 int Apply(svtkDataObject *dobj, DatasetFunction &func);
180 
181 /// Store ghost layer metadata in the mesh
182 SENSEI_EXPORT
183 int SetGhostLayerMetadata(svtkDataObject *mesh,
184  int nGhostCellLayers, int nGhostNodeLayers);
185 
186 /** Retreive ghost layer metadata from the mesh. returns non-zero if no such
187 * metadata is found.
188 */
189 SENSEI_EXPORT
190 int GetGhostLayerMetadata(svtkDataObject *mesh,
191  int &nGhostCellLayers, int &nGhostNodeLayers);
192 
193 /*** Get metadata, note that data set variant is not meant to be used on blocks
194  * of a multi-block
195  */
196 SENSEI_EXPORT
197 int GetMetadata(MPI_Comm comm, svtkDataSet *ds, MeshMetadataPtr);
198 SENSEI_EXPORT
199 int GetMetadata(MPI_Comm comm, svtkCompositeDataSet *cd, MeshMetadataPtr);
200 
201 /** Given a data object ensure that it is a composite data set If it already is,
202  * then the call is a no-op, if it is not then it is converted to a multiblock.
203  * The flag take determines if the smart pointer takes ownership or adds a
204  * reference.
205  */
206 SENSEI_EXPORT
207 svtkCompositeDataSetPtr AsCompositeData(MPI_Comm comm,
208  svtkDataObject *dobj, bool take = true);
209 
210 /// Return true if the mesh or block type is AMR
211 inline bool AMR(const MeshMetadataPtr &md)
212 {
213  return (md->MeshType == SVTK_OVERLAPPING_AMR) ||
214  (md->MeshType == SVTK_NON_OVERLAPPING_AMR);
215 }
216 
217 /// Return true if the mesh or block type is logically Cartesian
218 inline bool Structured(const MeshMetadataPtr &md)
219 {
220  return (md->BlockType == SVTK_STRUCTURED_GRID) ||
221  (md->MeshType == SVTK_STRUCTURED_GRID);
222 }
223 
224 /// Return true if the mesh or block type is polydata
225 inline bool Polydata(const MeshMetadataPtr &md)
226 {
227  return (md->BlockType == SVTK_POLY_DATA) || (md->MeshType == SVTK_POLY_DATA);
228 }
229 
230 /// Return true if the mesh or block type is unstructured
231 inline bool Unstructured(const MeshMetadataPtr &md)
232 {
233  return (md->BlockType == SVTK_UNSTRUCTURED_GRID) ||
234  (md->MeshType == SVTK_UNSTRUCTURED_GRID);
235 }
236 
237 /// Return true if the mesh or block type is stretched Cartesian
238 inline bool StretchedCartesian(const MeshMetadataPtr &md)
239 {
240  return (md->BlockType == SVTK_RECTILINEAR_GRID) ||
241  (md->MeshType == SVTK_RECTILINEAR_GRID);
242 }
243 
244 /// Return true if the mesh or block type is uniform Cartesian
245 inline bool UniformCartesian(const MeshMetadataPtr &md)
246 {
247  return (md->BlockType == SVTK_IMAGE_DATA) || (md->MeshType == SVTK_IMAGE_DATA)
248  || (md->BlockType == SVTK_UNIFORM_GRID) || (md->MeshType == SVTK_UNIFORM_GRID);
249 }
250 
251 /// Return true if the mesh or block type is logically Cartesian
252 inline bool LogicallyCartesian(const MeshMetadataPtr &md)
253 {
254  return Structured(md) || UniformCartesian(md) || StretchedCartesian(md);
255 }
256 
257 // rank 0 writes a dataset for visualizing the domain decomp
258 SENSEI_EXPORT
259 int WriteDomainDecomp(MPI_Comm comm, const sensei::MeshMetadataPtr &md,
260  const std::string fileName);
261 
262 /** Write an svtkImageData data set to disk in a VTK compatible binary format.
263  * This format can be read by ParaView and VisIt amoung others. This
264  * implementation does not depend on VTK iteself and can be used when VTK is
265  * not present in the build.
266  *
267  * @param[in] fn the name of the file to write
268  * @param[in] npx the number of points in the x-dimension
269  * @param[in] npy the number of points in the y-dimension
270  * @param[in] npz the number of points in the z-dimension
271  * @param[in] x0 the origin in the x-dimension
272  * @param[in] y0 the origin in the y-dimension
273  * @param[in] z0 the origin in the x-dimension
274  * @param[in] dx the grid spacing in the x-dimension
275  * @param[in] dy the grid spacing in the y-dimension
276  * @param[in] dz the grid spacing in the z-dimension
277  * @param[in] cellData a vector of cell centered arrays to write
278  * @param[in] pointData a vector of point centered arrays to write
279  * @returns 0 if successful
280  */
281 SENSEI_EXPORT
282 int WriteVTK(const char *fn, long npx, long npy, long npz,
283  double x0, double y0, double z0, double dx, double dy, double dz,
284  const std::vector<svtkDataArray*> &cellData,
285  const std::vector<svtkDataArray*> &pointData);
286 
287 
288 /** Packs data from a cell array into another cell array keeping track of
289  * where to insert into the output array. Use it to serialze verys, lines, polys
290  * strips form a polytdata into a single cell array for transport
291  */
292 template <typename SVTK_TT, typename ARRAY_TT = svtkAOSDataArrayTemplate<SVTK_TT>>
293 void PackCells(ARRAY_TT *coIn, ARRAY_TT *ccIn, ARRAY_TT *coOut, ARRAY_TT *ccOut,
294  size_t &coId, size_t &ccId)
295 {
296  // copy offsets
297  size_t nOffs = coIn->GetNumberOfTuples();
298  const SVTK_TT *pSrc = coIn->GetPointer(0);
299  SVTK_TT *pDest = coOut->GetPointer(0);
300 
301  for (size_t i = 0; i < nOffs; ++i)
302  pDest[coId + i] = pSrc[i];
303 
304  // update cell id
305  coId += nOffs;
306 
307  // copy connectivity
308  size_t nConn = ccIn->GetNumberOfTuples();
309  pSrc = ccIn->GetPointer(0);
310  pDest = ccOut->GetPointer(0);
311 
312  for (size_t i = 0; i < nConn; ++i)
313  pDest[ccId + i] = pSrc[i];
314 
315  // update conn id
316  ccId += nConn;
317 }
318 
319 /// deserializes a buffer made by SVTKUtils::PackCells.
320 template <typename SVTK_TT, typename ARRAY_TT = svtkAOSDataArrayTemplate<SVTK_TT>>
321 void UnpackCells(size_t nc, ARRAY_TT *coIn, ARRAY_TT *ccIn,
322  svtkCellArray *caOut, size_t &coId, size_t &ccId)
323 {
324  // offsets
325  size_t nCo = nc + 1;
326 
327  SVTK_TT *pCoIn = coIn->GetPointer(0);
328 
329  ARRAY_TT *coOut = ARRAY_TT::New();
330 
331  coOut->SetNumberOfTuples(nCo);
332  SVTK_TT *pCoOut = coOut->GetPointer(0);
333 
334  for (size_t i = 0; i < nCo; ++i)
335  pCoOut[i] = pCoIn[coId + i];
336 
337  coId += nCo;
338 
339  // connectivity
340  size_t nCc = nc ? pCoOut[nCo - 1] : 0;
341 
342  ARRAY_TT *ccOut = ARRAY_TT::New();
343  if (nCc)
344  {
345  SVTK_TT *pCcIn = ccIn->GetPointer(0);
346 
347  ccOut->SetNumberOfTuples(nCc);
348  SVTK_TT *pCcOut = ccOut->GetPointer(0);
349 
350  for (size_t i = 0; i < nCc; ++i)
351  pCcOut[i] = pCcIn[ccId + i];
352 
353  ccId += nCc;
354  }
355 
356  // package as cells
357  caOut->SetData(coOut, ccOut);
358 
359  coOut->Delete();
360  ccOut->Delete();
361 }
362 
363 /** Constructs VTK objects from SVTK objects enabling the use of VTK filters
364  * and ParaView Catalyst on SVTK data. The factory currently supports a
365  * limitted number of VTK objects but can be expanded as needed.
366  * See sensei::SVTKUtils::SVTKObjectFactory for conversions from VTK to SVTK.
367  */
369 {
370 public:
371  /** Construct a new VTK vtkDataArray from the passsed SVTK svtkDataArray.
372  * Returns a newly allocated instance of the corresponding VTK data array.
373  * Data is zero-copy transfered. A references to the passed SVTK
374  * svtkDataArray held by the newly cretaed VTK vtkDataArray ensuring
375  * propper life time. It is the callers responsibility to Delete the
376  * returned vtkDataArray instance when finished.
377  */
378  static vtkDataArray *New(svtkDataArray *daIn);
379 
380  /// overload for 64 bit cell arrays
381  static vtkTypeInt64Array *New(svtkTypeInt64Array *daIn);
382 
383  /// overload for 32 bit cell arrays
384  static vtkTypeInt32Array *New(svtkTypeInt32Array *daIn);
385 
386  /** Construct a new VTK vtkDataObject from the passsed SVTK svtkDataObject.
387  * Returns a newly allocated instance of the corresponding VTK object.
388  * Embedded svtkDataArray instances are zero-copy transfered. References to
389  * vtkDataArray transfered are held by the VTK object ensuring propper life
390  * time. It is the callers responsibility to Delete the returned
391  * vtkDataObject when finished. See the overloaded New methods for a list of
392  * implemented data objects.
393  */
394  static vtkDataObject *New(svtkDataObject *objIn);
395 
396  /// @copydoc New(svtkDataObject*)
397  static vtkDataSet *New(svtkDataSet *dsIn);
398 
399  /// @copydoc New(svtkDataObject*)
400  static vtkCellArray *New(svtkCellArray *caIn);
401 
402  /// @copydoc New(svtkDataObject*)
403  static vtkFieldData *New(svtkFieldData *fdIn);
404 
405  /// @copydoc New(svtkDataObject*)
406  static vtkPointData *New(svtkPointData *fdIn);
407 
408  /// @copydoc New(svtkDataObject*)
409  static vtkCellData *New(svtkCellData *fdIn);
410 
411  /// @copydoc New(svtkDataObject*)
412  static vtkPoints *New(svtkPoints *ptsIn);
413 
414  /// @copydoc New(svtkDataObject*)
415  static vtkImageData *New(svtkImageData *idIn);
416 
417  /// @copydoc New(svtkDataObject*)
418  static vtkUniformGrid *New(svtkUniformGrid *idIn);
419 
420  /// @copydoc New(svtkDataObject*)
421  static vtkRectilinearGrid *New(svtkRectilinearGrid *rgIn);
422 
423  /// @copydoc New(svtkDataObject*)
424  static vtkStructuredGrid *New(svtkStructuredGrid *sgIn);
425 
426  /// @copydoc New(svtkDataObject*)
427  static vtkPolyData *New(svtkPolyData *pdIn);
428 
429  /// @copydoc New(svtkDataObject*)
430  static vtkUnstructuredGrid *New(svtkUnstructuredGrid *ugIn);
431 
432  /// @copydoc New(svtkDataObject*)
433  static vtkMultiBlockDataSet *New(svtkMultiBlockDataSet *mbIn);
434 
435  /// @copydoc New(svtkDataObject*)
436  static vtkOverlappingAMR *New(svtkOverlappingAMR *amrIn);
437 };
438 
439 
440 
441 /** Constructs SVTK objects from VTK objects enabling the consumption of the
442  * output of VTK filters and ParaView Catalyst. The factory currently supports
443  * a limitted number of SVTK objects but can be expanded as needed.
444  * See sensei::SVTKUtils::VTKObjectFactory for conversions from SVTK to VTK.
445  */
447 {
448 public:
449  /** Construct a new SVTK svtkDataArray from the passsed VTK vtkDataArray.
450  * Returns a newly allocated instance of the corresponding SVTK data array.
451  * Data is zero-copy transfered. A references to the passed VTK
452  * vtkDataArray held by the newly cretaed SVTK svtkDataArray ensuring
453  * propper life time. It is the callers responsibility to Delete the
454  * returned svtkDataArray instance when finished.
455  */
456  static svtkDataArray *New(vtkDataArray *daIn);
457 
458  /// overload for 64 bit cell arrays
459  static svtkTypeInt64Array *New(vtkTypeInt64Array *daIn);
460 
461  /// overload for 32 bit cell arrays
462  static svtkTypeInt32Array *New(vtkTypeInt32Array *daIn);
463 
464  /** Construct a new SVTK svtkDataObject from the passsed VTK vtkDataObject.
465  * Returns a newly allocated instance of the corresponding SVTK object.
466  * Embedded vtkDataArray instances are zero-copy transfered. References to
467  * svtkDataArray transfered are held by the SVTK object ensuring propper life
468  * time. It is the callers responsibility to Delete the returned
469  * svtkDataObject when finished. See the overloaded New methods for a list of
470  * implemented data objects.
471  */
472  static svtkDataObject *New(vtkDataObject *objIn);
473 
474  /// @copydoc New(vtkDataObject*)
475  static svtkDataSet *New(vtkDataSet *dsIn);
476 
477  /// @copydoc New(vtkDataObject*)
478  static svtkCellArray *New(vtkCellArray *caIn);
479 
480  /// @copydoc New(vtkDataObject*)
481  static svtkFieldData *New(vtkFieldData *fdIn);
482 
483  /// @copydoc New(vtkDataObject*)
484  static svtkPointData *New(vtkPointData *fdIn);
485 
486  /// @copydoc New(vtkDataObject*)
487  static svtkCellData *New(vtkCellData *fdIn);
488 
489  /// @copydoc New(vtkDataObject*)
490  static svtkPoints *New(vtkPoints *ptsIn);
491 
492  /// @copydoc New(vtkDataObject*)
493  static svtkImageData *New(vtkImageData *idIn);
494 
495  /// @copydoc New(vtkDataObject*)
496  static svtkUniformGrid *New(vtkUniformGrid *idIn);
497 
498  /// @copydoc New(vtkDataObject*)
499  static svtkRectilinearGrid *New(vtkRectilinearGrid *rgIn);
500 
501  /// @copydoc New(vtkDataObject*)
502  static svtkStructuredGrid *New(vtkStructuredGrid *sgIn);
503 
504  /// @copydoc New(vtkDataObject*)
505  static svtkPolyData *New(vtkPolyData *pdIn);
506 
507  /// @copydoc New(vtkDataObject*)
508  static svtkUnstructuredGrid *New(vtkUnstructuredGrid *ugIn);
509 
510  /// @copydoc New(vtkDataObject*)
511  static svtkMultiBlockDataSet *New(vtkMultiBlockDataSet *mbIn);
512 
513  /// @copydoc New(vtkDataObject*)
514  static svtkOverlappingAMR *New(vtkOverlappingAMR *amrIn);
515 };
516 
517 }
518 }
519 
520 #endif
SENSEI_EXPORT const char * GetAttributesName(int association)
returns the name of the association, point, cell or field
bool Structured(const MeshMetadataPtr &md)
Return true if the mesh or block type is logically Cartesian.
Definition: SVTKUtils.h:218
SENSEI_EXPORT int GetGhostLayerMetadata(svtkDataObject *mesh, int &nGhostCellLayers, int &nGhostNodeLayers)
Retreive ghost layer metadata from the mesh.
bool StretchedCartesian(const MeshMetadataPtr &md)
Return true if the mesh or block type is stretched Cartesian.
Definition: SVTKUtils.h:238
bool AMR(const MeshMetadataPtr &md)
Return true if the mesh or block type is AMR.
Definition: SVTKUtils.h:211
static vtkDataArray * New(svtkDataArray *daIn)
Construct a new VTK vtkDataArray from the passsed SVTK svtkDataArray.
SENSEI_EXPORT unsigned int Size(int svtkt)
given a SVTK type enum returns the sizeof that type
Constructs VTK objects from SVTK objects enabling the use of VTK filters and ParaView Catalyst on SVT...
Definition: SVTKUtils.h:368
bool UniformCartesian(const MeshMetadataPtr &md)
Return true if the mesh or block type is uniform Cartesian.
Definition: SVTKUtils.h:245
SENSEI_EXPORT svtkCompositeDataSetPtr AsCompositeData(MPI_Comm comm, svtkDataObject *dobj, bool take=true)
Given a data object ensure that it is a composite data set If it already is, then the call is a no-op...
SENSEI_EXPORT int WriteVTK(const char *fn, long npx, long npy, long npz, double x0, double y0, double z0, double dx, double dy, double dz, const std::vector< svtkDataArray *> &cellData, const std::vector< svtkDataArray *> &pointData)
Write an svtkImageData data set to disk in a VTK compatible binary format.
void PackCells(ARRAY_TT *coIn, ARRAY_TT *ccIn, ARRAY_TT *coOut, ARRAY_TT *ccOut, size_t &coId, size_t &ccId)
Packs data from a cell array into another cell array keeping track of where to insert into the output...
Definition: SVTKUtils.h:293
void UnpackCells(size_t nc, ARRAY_TT *coIn, ARRAY_TT *ccIn, svtkCellArray *caOut, size_t &coId, size_t &ccId)
deserializes a buffer made by SVTKUtils::PackCells.
Definition: SVTKUtils.h:321
SENSEI_EXPORT svtkFieldData * GetAttributes(svtkDataSet *dobj, int association)
returns the container for the associations: svtkPointData, svtkCellData, or svtkFieldData ...
bool Unstructured(const MeshMetadataPtr &md)
Return true if the mesh or block type is unstructured.
Definition: SVTKUtils.h:231
Constructs SVTK objects from VTK objects enabling the consumption of the output of VTK filters and Pa...
Definition: SVTKUtils.h:446
SENSEI.
Definition: ADIOS2AnalysisAdaptor.h:27
SENSEI_EXPORT int IsLegacyDataObject(int code)
given a SVTK data object enum returns true if it a legacy object
SENSEI_EXPORT int Apply(svtkDataObject *input, svtkDataObject *output, BinaryDatasetFunction &func)
Applies the function to leaves of the structurally equivalent input and output data objects...
bool Polydata(const MeshMetadataPtr &md)
Return true if the mesh or block type is polydata.
Definition: SVTKUtils.h:225
std::function< int(svtkDataSet *, svtkDataSet *)> BinaryDatasetFunction
A callback that processes input and output datasets.
Definition: SVTKUtils.h:161
std::function< int(svtkDataSet *)> DatasetFunction
A callback that processes input and output datasets return 0 for success, > zero to stop without erro...
Definition: SVTKUtils.h:173
SENSEI_EXPORT svtkDataObject * NewDataObject(int code)
givne a SVTK data object enum constructs an instance
SENSEI_EXPORT int SetGhostLayerMetadata(svtkDataObject *mesh, int nGhostCellLayers, int nGhostNodeLayers)
Store ghost layer metadata in the mesh.
bool LogicallyCartesian(const MeshMetadataPtr &md)
Return true if the mesh or block type is logically Cartesian.
Definition: SVTKUtils.h:252
SENSEI_EXPORT int GetAssociation(std::string assocStr, int &assoc)
returns the enum value given an association name.
SVTK_TT * GetPointer(svtkDataArray *da)
given a svtkDataArray get a pointer to underlying data this handles access from SVTK&#39;s AOS and SOA la...
Definition: SVTKUtils.h:108