File: vtkHyperTreeGridToUnstructuredGrid.cxx

package info (click to toggle)
vtk9 9.5.2%2Bdfsg3-4
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 205,916 kB
  • sloc: cpp: 2,336,565; ansic: 327,116; python: 111,200; yacc: 4,104; java: 3,977; sh: 3,032; xml: 2,771; perl: 2,189; lex: 1,787; makefile: 178; javascript: 165; objc: 153; tcl: 59
file content (289 lines) | stat: -rw-r--r-- 8,233 bytes parent folder | download | duplicates (7)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#include "vtkHyperTreeGridToUnstructuredGrid.h"

#include "vtkBitArray.h"
#include "vtkCellArray.h"
#include "vtkCellData.h"
#include "vtkHyperTreeGrid.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkUnstructuredGrid.h"

#include "vtkHyperTreeGridNonOrientedGeometryCursor.h"

VTK_ABI_NAMESPACE_BEGIN
vtkStandardNewMacro(vtkHyperTreeGridToUnstructuredGrid);

//------------------------------------------------------------------------------
vtkHyperTreeGridToUnstructuredGrid::vtkHyperTreeGridToUnstructuredGrid()
  : Points(nullptr)
  , Cells(nullptr)
  , Dimension(0)
  , Orientation(0)
  , Axes(nullptr)
  , AddOriginalIds(true)
  , OriginalIds(nullptr)
{
}

//------------------------------------------------------------------------------
// The class members are only used during process and are destroyed once the
// process is finished to reduce stack size during recursive calls.
vtkHyperTreeGridToUnstructuredGrid::~vtkHyperTreeGridToUnstructuredGrid() = default;

//------------------------------------------------------------------------------
void vtkHyperTreeGridToUnstructuredGrid::PrintSelf(ostream& os, vtkIndent indent)
{
  this->Superclass::PrintSelf(os, indent);
}

//------------------------------------------------------------------------------
int vtkHyperTreeGridToUnstructuredGrid::FillOutputPortInformation(int, vtkInformation* info)
{
  info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
  return 1;
}

//------------------------------------------------------------------------------
int vtkHyperTreeGridToUnstructuredGrid::ProcessTrees(
  vtkHyperTreeGrid* input, vtkDataObject* outputDO)
{
  // Downcast output data object to hyper tree grid
  vtkUnstructuredGrid* output = vtkUnstructuredGrid::SafeDownCast(outputDO);
  if (!output)
  {
    vtkErrorMacro("Incorrect type of output: " << outputDO->GetClassName());
    return 0;
  }

  // Set instance variables needed for this conversion
  this->Points = vtkPoints::New();
  this->Cells = vtkCellArray::New();
  this->Dimension = input->GetDimension();
  this->Orientation = input->GetOrientation();
  this->Axes = input->GetAxes();

  // Initialize output cell data
  this->InData = input->GetCellData();
  this->OutData = output->GetCellData();
  this->OutData->CopyAllocate(this->InData);

  if (this->AddOriginalIds)
  {
    this->OriginalIds = vtkIdTypeArray::New();
    this->OriginalIds->SetName("OriginalIds");
    this->OriginalIds->SetNumberOfComponents(1);
    this->OriginalIds->SetNumberOfTuples(input->GetNumberOfLeaves());
  }

  // Iterate over all hyper trees
  vtkIdType index;
  vtkHyperTreeGrid::vtkHyperTreeGridIterator it;
  input->InitializeTreeIterator(it);
  vtkNew<vtkHyperTreeGridNonOrientedGeometryCursor> cursor;
  while (it.GetNextTree(index))
  {
    if (this->CheckAbort())
    {
      break;
    }
    // Initialize new geometric cursor at root of current tree
    input->InitializeNonOrientedGeometryCursor(cursor, index);

    // Convert hyper tree into unstructured mesh recursively
    this->RecursivelyProcessTree(cursor);
  } // it

  // Set output geometry and topology
  output->SetPoints(this->Points);
  switch (this->Dimension)
  {
    case 1:
      // 1D cells are lines
      output->SetCells(VTK_LINE, this->Cells);
      break;
    case 2:
      // 2D cells are quadrilaterals
      output->SetCells(VTK_PIXEL, this->Cells);
      break;
    case 3:
      // 3D cells are voxels (i.e. hexahedra with indexing order equal to that of cursors)
      output->SetCells(VTK_VOXEL, this->Cells);
      break;
    default:
      break;
  } // switch ( this->Dimension )

  if (this->AddOriginalIds)
  {
    this->OutData->AddArray(this->OriginalIds);
    this->OriginalIds->FastDelete();
    this->OriginalIds = nullptr;
  }

  this->Points->FastDelete();
  this->Cells->FastDelete();
  this->Points = nullptr;
  this->Cells = nullptr;

  return 1;
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToUnstructuredGrid::RecursivelyProcessTree(
  vtkHyperTreeGridNonOrientedGeometryCursor* cursor)
{
  // If leaf is masked, skip it
  if (cursor->IsMasked())
  {
    return;
  }

  // Create unstructured output if cursor is at leaf
  if (cursor->IsLeaf())
  {
    // Cursor is at leaf, retrieve its global index
    vtkIdType id = cursor->GetGlobalNodeIndex();

    // Create cell
    this->AddCell(id, cursor->GetOrigin(), cursor->GetSize());
  } // if ( cursor->IsLeaf() )
  else
  {
    // Cursor is not at leaf, recurse to all children
    int numChildren = cursor->GetNumberOfChildren();
    for (int ichild = 0; ichild < numChildren; ++ichild)
    {
      if (this->CheckAbort())
      {
        break;
      }
      cursor->ToChild(ichild);
      // Recurse
      this->RecursivelyProcessTree(cursor);
      cursor->ToParent();
    } // child
  }   // else
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToUnstructuredGrid::AddCell(vtkIdType inId, double* origin, double* size)
{
  // Storage for point coordinates
  double pt[] = { 0., 0., 0. };

  // Storage for cell vertex IDs
  vtkIdType ids[8];

  // Storage for cell ID
  vtkIdType outId;

  // First cell vertex is always at origin of cursor
  // Add vertex #0 : (0,0)
  memcpy(pt, origin, 3 * sizeof(double));
  ids[0] = this->Points->InsertNextPoint(pt);

  // Create remaining 2^d - 1 vertices depending on dimension
  switch (this->Dimension)
  {
    case 1:
    {
      assert("pre: internal" && this->Orientation == this->Axes[0]);

      // In 1D there is only one other vertex
      pt[0] = origin[this->Orientation] + size[this->Orientation];
      ids[1] = this->Points->InsertNextPoint(pt);

      // Insert next line
      outId = this->Cells->InsertNextCell(2, ids);
      break;
    }
    case 2:
    {
      unsigned int axis1 = this->Axes[0];
      unsigned int axis2 = this->Axes[1];

      // Add vertex #1 : (1,0)
      pt[axis1] = origin[axis1] + size[axis1];
      pt[axis2] = origin[axis2];
      ids[1] = this->Points->InsertNextPoint(pt);

      // Add vertex #2 : (0,1)
      pt[axis1] = origin[axis1];
      pt[axis2] = origin[axis2] + size[axis2];
      ids[2] = this->Points->InsertNextPoint(pt);

      // Add vertex #3 : (1,1)
      pt[axis1] = origin[axis1] + size[axis1];
      pt[axis2] = origin[axis2] + size[axis2];
      ids[3] = this->Points->InsertNextPoint(pt);

      // Insert next quadrangle
      outId = this->Cells->InsertNextCell(4, ids);
      break;
    }
    case 3:
    {
      // z=0 plane
      pt[2] = origin[2];

      // Add vertex #1 : (1,0,0)
      pt[0] = origin[0] + size[0];
      pt[1] = origin[1];
      ids[1] = this->Points->InsertNextPoint(pt);

      // Add vertex #2 : (0,1,0)
      pt[0] = origin[0];
      pt[1] = origin[1] + size[1];
      ids[2] = this->Points->InsertNextPoint(pt);

      // Add vertex #3 : (1,1,0)
      pt[0] = origin[0] + size[0];
      pt[1] = origin[1] + size[1];
      ids[3] = this->Points->InsertNextPoint(pt);

      // z=1 plane
      pt[2] = origin[2] + size[2];

      // Add vertex #4 : (0,0,1)
      pt[0] = origin[0];
      pt[1] = origin[1];
      ids[4] = this->Points->InsertNextPoint(pt);

      // Add vertex #5 : (1,0,1)
      pt[0] = origin[0] + size[0];
      pt[1] = origin[1];
      ids[5] = this->Points->InsertNextPoint(pt);

      // Add vertex #6 : (0,1,1)
      pt[0] = origin[0];
      pt[1] = origin[1] + size[1];
      ids[6] = this->Points->InsertNextPoint(pt);

      // Add vertex #7 : (1,1,1)
      pt[0] = origin[0] + size[0];
      pt[1] = origin[1] + size[1];
      ids[7] = this->Points->InsertNextPoint(pt);

      // Insert next voxel
      outId = this->Cells->InsertNextCell(8, ids);
      break;
    }
    default:
    {
      return;
    }
  } // switch ( this->Dimension )

  // Copy output data from input
  this->OutData->CopyData(this->InData, inId, outId);

  // And the global id if needed
  if (this->AddOriginalIds)
  {
    this->OriginalIds->SetTuple1(outId, inId);
  }
}
VTK_ABI_NAMESPACE_END