VTK  9.0.1
vtkStructuredData.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkStructuredData.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 =========================================================================*/
29 #ifndef vtkStructuredData_h
30 #define vtkStructuredData_h
31 
32 #include "vtkCommonDataModelModule.h" // For export macro
33 #include "vtkObject.h"
34 
35 class vtkIdList;
36 
37 #define VTK_UNCHANGED 0
38 #define VTK_SINGLE_POINT 1
39 #define VTK_X_LINE 2
40 #define VTK_Y_LINE 3
41 #define VTK_Z_LINE 4
42 #define VTK_XY_PLANE 5
43 #define VTK_YZ_PLANE 6
44 #define VTK_XZ_PLANE 7
45 #define VTK_XYZ_GRID 8
46 #define VTK_EMPTY 9
47 
48 class VTKCOMMONDATAMODEL_EXPORT vtkStructuredData : public vtkObject
49 {
50 public:
51  vtkTypeMacro(vtkStructuredData, vtkObject);
52 
54 
61  static int SetDimensions(int inDim[3], int dim[3]);
62  static int SetExtent(int inExt[6], int ext[6]);
64 
66 
70  static int GetDataDescription(int dims[3]);
71  static int GetDataDescriptionFromExtent(int ext[6]);
73 
75 
78  static int GetDataDimension(int dataDescription);
79  static int GetDataDimension(int ext[6]);
81 
87  static vtkIdType GetNumberOfPoints(const int ext[6], int dataDescription = VTK_EMPTY);
88 
94  static vtkIdType GetNumberOfCells(const int ext[6], int dataDescription = VTK_EMPTY);
95 
101  static void GetCellExtentFromPointExtent(
102  const int pntExtent[6], int cellExtent[6], int dataDescription = VTK_EMPTY);
103 
108  static void GetDimensionsFromExtent(
109  const int ext[6], int dims[3], int dataDescription = VTK_EMPTY);
110 
117  static void GetCellDimensionsFromExtent(
118  const int ext[6], int celldims[3], int dataDescription = VTK_EMPTY);
119 
125  static void GetCellDimensionsFromPointDimensions(const int pntdims[3], int cellDims[3]);
126 
133  static void GetLocalStructuredCoordinates(
134  const int ijk[3], const int ext[6], int lijk[3], int dataDescription = VTK_EMPTY);
135 
141  static void GetGlobalStructuredCoordinates(
142  const int lijk[3], const int ext[6], int ijk[3], int dataDescription = VTK_EMPTY);
143 
147  static void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds, int dataDescription, int dim[3]);
148 
152  static void GetPointCells(vtkIdType ptId, vtkIdList* cellIds, int dim[3]);
153 
158  static void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int dim[3]);
159  static void GetCellNeighbors(
160  vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds, int dim[3], int seedLoc[3]);
161 
167  static vtkIdType ComputePointIdForExtent(
168  const int extent[6], const int ijk[3], int dataDescription = VTK_EMPTY);
169 
175  static vtkIdType ComputeCellIdForExtent(
176  const int extent[6], const int ijk[3], int dataDescription = VTK_EMPTY);
177 
184  static vtkIdType ComputePointId(
185  const int dim[3], const int ijk[3], int dataDescription = VTK_EMPTY);
186 
193  static vtkIdType ComputeCellId(
194  const int dim[3], const int ijk[3], int dataDescription = VTK_EMPTY);
195 
202  static void ComputeCellStructuredCoordsForExtent(
203  const vtkIdType cellIdx, const int ext[6], int ijk[3], int dataDescription = VTK_EMPTY);
204 
210  static void ComputeCellStructuredCoords(
211  const vtkIdType cellId, const int dim[3], int ijk[3], int dataDescription = VTK_EMPTY);
212 
218  static void ComputePointStructuredCoordsForExtent(
219  const vtkIdType ptId, const int ext[6], int ijk[3], int dataDescription = VTK_EMPTY);
220 
226  static void ComputePointStructuredCoords(
227  const vtkIdType ptId, const int dim[3], int ijk[3], int dataDescription = VTK_EMPTY);
228 
229 protected:
231  ~vtkStructuredData() override {}
232 
240  static vtkIdType GetLinearIndex(const int i, const int j, const int k, const int N1, const int N2)
241  {
242  return ((static_cast<vtkIdType>(k) * N2 + j) * N1 + i);
243  }
244 
246 
253  const vtkIdType idx, const int N1, const int N2, int& i, int& j, int& k)
254  {
255  vtkIdType N12 = N1 * N2;
256  k = static_cast<int>(idx / N12);
257  j = static_cast<int>((idx - k * N12) / N1);
258  i = static_cast<int>(idx - k * N12 - j * N1);
259  }
261 
262  // Want to avoid importing <algorithm> in the header...
263  template <typename T>
264  static T Max(const T& a, const T& b)
265  {
266  return (a > b) ? a : b;
267  }
268 
269 private:
270  vtkStructuredData(const vtkStructuredData&) = delete;
271  void operator=(const vtkStructuredData&) = delete;
272 };
273 
274 //------------------------------------------------------------------------------
275 inline void vtkStructuredData::GetCellDimensionsFromExtent(const int ext[6], int celldims[3], int)
276 {
277  celldims[0] = vtkStructuredData::Max(ext[1] - ext[0], 0);
278  celldims[1] = vtkStructuredData::Max(ext[3] - ext[2], 0);
279  celldims[2] = vtkStructuredData::Max(ext[5] - ext[4], 0);
280 }
281 
282 //------------------------------------------------------------------------------
283 inline vtkIdType vtkStructuredData::ComputePointId(const int dims[3], const int ijk[3], int)
284 {
285  return vtkStructuredData::GetLinearIndex(ijk[0], ijk[1], ijk[2], dims[0], dims[1]);
286 }
287 
288 //------------------------------------------------------------------------------
289 inline vtkIdType vtkStructuredData::ComputeCellId(const int dims[3], const int ijk[3], int)
290 {
291  return vtkStructuredData::GetLinearIndex(ijk[0], ijk[1], ijk[2],
292  vtkStructuredData::Max(dims[0] - 1, 1), vtkStructuredData::Max(dims[1] - 1, 1));
293 }
294 
295 //------------------------------------------------------------------------------
296 inline vtkIdType vtkStructuredData::GetNumberOfPoints(const int ext[6], int)
297 {
298  return static_cast<vtkIdType>(ext[1] - ext[0] + 1) * static_cast<vtkIdType>(ext[3] - ext[2] + 1) *
299  static_cast<vtkIdType>(ext[5] - ext[4] + 1);
300 }
301 
302 //------------------------------------------------------------------------------
303 inline vtkIdType vtkStructuredData::GetNumberOfCells(const int ext[6], int)
304 {
305  int cellDims[3];
307 
308  // Replace 0's with 1's so we can just multiply them regardless of cell type.
309  cellDims[0] = vtkStructuredData::Max(cellDims[0], 1);
310  cellDims[1] = vtkStructuredData::Max(cellDims[1], 1);
311  cellDims[2] = vtkStructuredData::Max(cellDims[2], 1);
312 
313  // Note, when we compute the result below, we statically cast to vtkIdType to
314  // ensure the compiler will generate a 32x32=64 instruction.
315  return static_cast<vtkIdType>(cellDims[0]) * static_cast<vtkIdType>(cellDims[1]) *
316  static_cast<vtkIdType>(cellDims[2]);
317 }
318 
319 //------------------------------------------------------------------------------
321  const int nodeExtent[6], int cellExtent[6], int)
322 {
323  cellExtent[0] = nodeExtent[0];
324  cellExtent[2] = nodeExtent[2];
325  cellExtent[4] = nodeExtent[4];
326 
327  cellExtent[1] = vtkStructuredData::Max(nodeExtent[0], nodeExtent[1] - 1);
328  cellExtent[3] = vtkStructuredData::Max(nodeExtent[2], nodeExtent[3] - 1);
329  cellExtent[5] = vtkStructuredData::Max(nodeExtent[4], nodeExtent[5] - 1);
330 }
331 
332 //------------------------------------------------------------------------------
333 inline void vtkStructuredData::GetDimensionsFromExtent(const int ext[6], int dims[3], int)
334 {
335  dims[0] = ext[1] - ext[0] + 1;
336  dims[1] = ext[3] - ext[2] + 1;
337  dims[2] = ext[5] - ext[4] + 1;
338 }
339 
340 //------------------------------------------------------------------------------
342  const int nodeDims[3], int cellDims[3])
343 {
344  cellDims[0] = vtkStructuredData::Max(nodeDims[0] - 1, 0);
345  cellDims[1] = vtkStructuredData::Max(nodeDims[1] - 1, 0);
346  cellDims[2] = vtkStructuredData::Max(nodeDims[2] - 1, 0);
347 }
348 
349 //------------------------------------------------------------------------------
351  const int ijk[3], const int ext[6], int lijk[3], int)
352 {
353  lijk[0] = ijk[0] - ext[0];
354  lijk[1] = ijk[1] - ext[2];
355  lijk[2] = ijk[2] - ext[4];
356 }
357 
358 //------------------------------------------------------------------------------
360  const int lijk[3], const int ext[6], int ijk[3], int)
361 {
362  ijk[0] = ext[0] + lijk[0];
363  ijk[1] = ext[2] + lijk[1];
364  ijk[2] = ext[4] + lijk[2];
365 }
366 
367 //------------------------------------------------------------------------------
369  const int extent[6], const int ijk[3], int)
370 {
371  int dims[3];
373 
374  int lijk[3];
376 
377  return vtkStructuredData::ComputePointId(dims, lijk);
378 }
379 
380 //------------------------------------------------------------------------------
382  const int extent[6], const int ijk[3], int)
383 {
384  int nodeDims[3];
386 
387  int lijk[3];
389 
390  return vtkStructuredData::ComputeCellId(nodeDims, lijk);
391 }
392 
393 //------------------------------------------------------------------------------
395  const vtkIdType cellId, const int dims[3], int ijk[3], int)
396 {
398  cellId, dims[0] - 1, dims[1] - 1, ijk[0], ijk[1], ijk[2]);
399 }
400 
401 //------------------------------------------------------------------------------
403  const vtkIdType cellIdx, const int ext[6], int ijk[3], int)
404 {
405  int nodeDims[3];
407 
408  int lijk[3];
409  vtkStructuredData::ComputeCellStructuredCoords(cellIdx, nodeDims, lijk);
410 
412 }
413 
414 //------------------------------------------------------------------------------
416  const vtkIdType ptId, const int dim[3], int ijk[3], int)
417 {
418  vtkStructuredData::GetStructuredCoordinates(ptId, dim[0], dim[1], ijk[0], ijk[1], ijk[2]);
419 }
420 
421 //------------------------------------------------------------------------------
423  const vtkIdType ptId, const int ext[6], int ijk[3], int)
424 {
425  int nodeDims[3];
427 
428  int lijk[3];
430 
432 }
433 
434 #endif
435 
436 // VTK-HeaderTest-Exclude: vtkStructuredData.h
vtkStructuredData::GetLocalStructuredCoordinates
static void GetLocalStructuredCoordinates(const int ijk[3], const int ext[6], int lijk[3], int dataDescription=VTK_EMPTY)
Given the global structured coordinates for a point or cell, ijk, w.r.t.
Definition: vtkStructuredData.h:350
vtkStructuredData::GetCellExtentFromPointExtent
static void GetCellExtentFromPointExtent(const int pntExtent[6], int cellExtent[6], int dataDescription=VTK_EMPTY)
Given the point extent of a grid, this method computes the corresponding cell extent for the grid.
Definition: vtkStructuredData.h:320
vtkIdType
int vtkIdType
Definition: vtkType.h:338
vtkStructuredData::GetCellDimensionsFromExtent
static void GetCellDimensionsFromExtent(const int ext[6], int celldims[3], int dataDescription=VTK_EMPTY)
Returns the cell dimensions, i.e., the number of cells along the i,j,k for the grid with the given gr...
Definition: vtkStructuredData.h:275
vtkStructuredData::ComputeCellStructuredCoords
static void ComputeCellStructuredCoords(const vtkIdType cellId, const int dim[3], int ijk[3], int dataDescription=VTK_EMPTY)
Given a cellId and grid dimensions 'dim', get the structured coordinates (i-j-k).
Definition: vtkStructuredData.h:394
vtkStructuredData::ComputePointStructuredCoordsForExtent
static void ComputePointStructuredCoordsForExtent(const vtkIdType ptId, const int ext[6], int ijk[3], int dataDescription=VTK_EMPTY)
Given a pointId and the grid extent ext, get the structured coordinates (i-j-k).
Definition: vtkStructuredData.h:422
vtkObject
abstract base class for most VTK objects
Definition: vtkObject.h:62
vtkStructuredData::GetDimensionsFromExtent
static void GetDimensionsFromExtent(const int ext[6], int dims[3], int dataDescription=VTK_EMPTY)
Computes the structured grid dimensions based on the given extent.
Definition: vtkStructuredData.h:333
vtkStructuredData::GetCellDimensionsFromPointDimensions
static void GetCellDimensionsFromPointDimensions(const int pntdims[3], int cellDims[3])
Given the dimensions of the grid, in pntdims, this method returns the corresponding cell dimensions f...
Definition: vtkStructuredData.h:341
vtkStructuredData::GetLinearIndex
static vtkIdType GetLinearIndex(const int i, const int j, const int k, const int N1, const int N2)
Computes the linear index for the given i-j-k structured of a grid with of N1 and N2 dimensions along...
Definition: vtkStructuredData.h:240
vtkStructuredData::GetNumberOfCells
static vtkIdType GetNumberOfCells(const int ext[6], int dataDescription=VTK_EMPTY)
Given the grid extent, this method returns the total number of cells within the extent.
Definition: vtkStructuredData.h:303
vtkStructuredData::GetGlobalStructuredCoordinates
static void GetGlobalStructuredCoordinates(const int lijk[3], const int ext[6], int ijk[3], int dataDescription=VTK_EMPTY)
Given local structured coordinates, and the corresponding global sub-grid extent, this method compute...
Definition: vtkStructuredData.h:359
vtkStructuredData::vtkStructuredData
vtkStructuredData()
Definition: vtkStructuredData.h:230
vtkStructuredData::~vtkStructuredData
~vtkStructuredData() override
Definition: vtkStructuredData.h:231
vtkStructuredData::Max
static T Max(const T &a, const T &b)
Definition: vtkStructuredData.h:264
vtkStructuredData::ComputePointId
static vtkIdType ComputePointId(const int dim[3], const int ijk[3], int dataDescription=VTK_EMPTY)
Given a location in structured coordinates (i-j-k), and the dimensions of the structured dataset,...
Definition: vtkStructuredData.h:283
vtkStructuredData
Singleton class for topologically regular data.
Definition: vtkStructuredData.h:48
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkStructuredData::GetNumberOfPoints
static vtkIdType GetNumberOfPoints(const int ext[6], int dataDescription=VTK_EMPTY)
Given the grid extent, this method returns the total number of points within the extent.
Definition: vtkStructuredData.h:296
vtkObject.h
vtkStructuredData::ComputePointStructuredCoords
static void ComputePointStructuredCoords(const vtkIdType ptId, const int dim[3], int ijk[3], int dataDescription=VTK_EMPTY)
Given a pointId and grid dimensions 'dim', get the structured coordinates (i-j-k).
Definition: vtkStructuredData.h:415
vtkStructuredData::ComputeCellStructuredCoordsForExtent
static void ComputeCellStructuredCoordsForExtent(const vtkIdType cellIdx, const int ext[6], int ijk[3], int dataDescription=VTK_EMPTY)
Given the global grid extent and the linear index of a cell within the grid extent,...
Definition: vtkStructuredData.h:402
vtkX3D::extent
@ extent
Definition: vtkX3D.h:351
vtkStructuredData::ComputeCellId
static vtkIdType ComputeCellId(const int dim[3], const int ijk[3], int dataDescription=VTK_EMPTY)
Given a location in structured coordinates (i-j-k), and the dimensions of the structured dataset,...
Definition: vtkStructuredData.h:289
vtkStructuredData::GetStructuredCoordinates
static void GetStructuredCoordinates(const vtkIdType idx, const int N1, const int N2, int &i, int &j, int &k)
Returns the structured coordinates (i,j,k) for the given linear index of a grid with N1 and N2 dimens...
Definition: vtkStructuredData.h:252
vtkStructuredData::ComputeCellIdForExtent
static vtkIdType ComputeCellIdForExtent(const int extent[6], const int ijk[3], int dataDescription=VTK_EMPTY)
Given a location in structured coordinates (i-j-k), and the extent of the structured dataset,...
Definition: vtkStructuredData.h:381
VTK_EMPTY
#define VTK_EMPTY
Definition: vtkStructuredData.h:46
vtkStructuredData::ComputePointIdForExtent
static vtkIdType ComputePointIdForExtent(const int extent[6], const int ijk[3], int dataDescription=VTK_EMPTY)
Given a location in structured coordinates (i-j-k), and the extent of the structured dataset,...
Definition: vtkStructuredData.h:368