File: vtkHyperTreeGridToDualGrid.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 (1310 lines) | stat: -rw-r--r-- 44,211 bytes parent folder | download | duplicates (6)
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
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#include "vtkHyperTreeGridToDualGrid.h"

#include "vtkBitArray.h"
#include "vtkCellData.h"
#include "vtkCellType.h"
#include "vtkHyperTree.h"
#include "vtkHyperTreeGrid.h"
#include "vtkHyperTreeGridNonOrientedMooreSuperCursor.h"
#include "vtkHyperTreeGridOrientedGeometryCursor.h"
#include "vtkIdTypeArray.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkPointData.h"
#include "vtkPoints.h"
#include "vtkUnstructuredGrid.h"

VTK_ABI_NAMESPACE_BEGIN
vtkStandardNewMacro(vtkHyperTreeGridToDualGrid);

static const unsigned int CornerNeighborCursorsTable3D0[8] = { 0, 1, 3, 4, 9, 10, 12, 13 };
static const unsigned int CornerNeighborCursorsTable3D1[8] = { 1, 2, 4, 5, 10, 11, 13, 14 };
static const unsigned int CornerNeighborCursorsTable3D2[8] = { 3, 4, 6, 7, 12, 13, 15, 16 };
static const unsigned int CornerNeighborCursorsTable3D3[8] = { 4, 5, 7, 8, 13, 14, 16, 17 };
static const unsigned int CornerNeighborCursorsTable3D4[8] = { 9, 10, 12, 13, 18, 19, 21, 22 };
static const unsigned int CornerNeighborCursorsTable3D5[8] = { 10, 11, 13, 14, 19, 20, 22, 23 };
static const unsigned int CornerNeighborCursorsTable3D6[8] = { 12, 13, 15, 16, 21, 22, 24, 25 };
static const unsigned int CornerNeighborCursorsTable3D7[8] = { 13, 14, 16, 17, 22, 23, 25, 26 };
static const unsigned int* CornerNeighborCursorsTable3D[8] = {
  CornerNeighborCursorsTable3D0,
  CornerNeighborCursorsTable3D1,
  CornerNeighborCursorsTable3D2,
  CornerNeighborCursorsTable3D3,
  CornerNeighborCursorsTable3D4,
  CornerNeighborCursorsTable3D5,
  CornerNeighborCursorsTable3D6,
  CornerNeighborCursorsTable3D7,
};

//------------------------------------------------------------------------------
vtkHyperTreeGridToDualGrid::vtkHyperTreeGridToDualGrid()
{
  // Dual grid corners (primal grid leaf centers)
  this->Points = nullptr;
  this->Connectivity = nullptr;
}

//------------------------------------------------------------------------------
vtkHyperTreeGridToDualGrid::~vtkHyperTreeGridToDualGrid()
{
  if (this->Points)
  {
    this->Points->Delete();
  }

  if (this->Connectivity)
  {
    this->Connectivity->Delete();
  }
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::PrintSelf(ostream& os, vtkIndent indent)
{
  this->Superclass::PrintSelf(os, indent);
  os << indent << "Points: " << this->Points << endl;
  os << indent << "Connectivity: " << this->Connectivity << endl;
}

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

//------------------------------------------------------------------------------
int vtkHyperTreeGridToDualGrid::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;
  }

  // Check if we can break out early
  // TODO: this isn't right, we need to trust the pipeline tell us when
  // we need to reexecute, or upstream grid will change and output of this filter will be stale.
  if (this->Points)
  {
    return 1;
  }

  // TODO: we shouldn't be using persistent internal Points and Connectivity, just populate the
  // output.
  // Create arrays needed by dual mesh
  this->Points = vtkPoints::New();
  this->Connectivity = vtkIdTypeArray::New();

  // Primal cell centers are dual points
  // TODO: We cannot only dimension the point array as the number of vertices,
  // especially if we want to keep a 1:1 mapping between HTG nodes and the dual mesh.
  // If we define a specific GlobalIndex or IndexStart, it would be too small,
  // because GetGlobalIndex returns a value greater than the number of cells.
  this->Points->SetNumberOfPoints(input->GetGlobalNodeIndexMax() + 1);

  // TODO: find out why we get some uninitialized point coords instead
  this->Points->GetData()->Fill(0.0);

  int numVerts = 1 << input->GetDimension();
  this->Connectivity->SetNumberOfComponents(numVerts);

  // Initialize grid depth
  unsigned int gridDepth = 0;

  // Compute and assign scales of all tree roots
  // double scale[3] = { 1., 1., 1. };

  // Check whether coordinate arrays match grid size
  // If coordinates array are complete, compute all tree scales
  // SEB: int* dims = input->GetDimensions();
  // SEB: if ( dims[0] == input->GetXCoordinates()->GetNumberOfTuples()
  // SEB:      && dims[1] == input->GetYCoordinates()->GetNumberOfTuples()
  // SEB:      && dims[2] == input->GetZCoordinates()->GetNumberOfTuples() )
  if (static_cast<int>(input->GetDimensions()[0]) ==
      input->GetXCoordinates()->GetNumberOfTuples() &&
    static_cast<int>(input->GetDimensions()[1]) == input->GetYCoordinates()->GetNumberOfTuples() &&
    static_cast<int>(input->GetDimensions()[2]) == input->GetZCoordinates()->GetNumberOfTuples())
  {
    gridDepth = input->GetNumberOfLevels();
  }

  // Compute and store reduction factors for speed
  double factor = 1.;
  for (unsigned int p = 0; p < gridDepth; ++p)
  {
    this->ReductionFactors[p] = .5 * factor;
    factor /= input->GetBranchFactor();
  } // p

  // Retrieve material mask
  vtkBitArray* mask = input->HasMask() ? input->GetMask() : nullptr;

  // Iterate over all hyper trees
  vtkIdType index;
  vtkHyperTreeGrid::vtkHyperTreeGridIterator it;
  input->InitializeTreeIterator(it);
  vtkNew<vtkHyperTreeGridNonOrientedMooreSuperCursor> cursor;
  while (it.GetNextTree(index))
  {
    if (this->CheckAbort())
    {
      break;
    }
    // Initialize new Moore cursor at root of current tree
    input->InitializeNonOrientedMooreSuperCursor(cursor, index);
    // Convert hyper tree into unstructured mesh recursively
    if (mask)
    {
      this->TraverseDualRecursively(cursor, mask, input);
    }
    else
    {
      this->TraverseDualRecursively(cursor, input);
    }
  } // it

  // Adjust dual points as needed to fit the primal boundary
  for (unsigned int d = 0; d < input->GetDimension(); ++d)
  {
    // Iterate over all adjustments for current dimension
    for (std::map<vtkIdType, double>::const_iterator _it = this->PointShifts[d].begin();
         _it != this->PointShifts[d].end(); ++_it)
    {
      double pt[3];

      assert(_it->first < input->GetNumberOfCells());
      this->Points->GetPoint(_it->first, pt);

      pt[d] += _it->second;

      assert(_it->first < input->GetNumberOfCells());
      this->Points->SetPoint(_it->first, pt);
    } // it
    this->PointShifts[d].clear();
  } // d
  this->PointShifted.clear();

  // now populate my output from the mesh internals made above
  output->SetPoints(this->Points);
  output->GetPointData()->ShallowCopy(input->GetCellData());

  int numPts = 1 << input->GetDimension();
  int type = VTK_VOXEL;
  if (numPts == 4)
  {
    type = VTK_PIXEL;
  }
  if (numPts == 2)
  {
    type = VTK_LINE;
  }
  output->Allocate();
  int numCells = this->Connectivity->GetNumberOfTuples();
  for (int cellIdx = 0; cellIdx < numCells; cellIdx++)
  {
    vtkIdType* ptr = this->Connectivity->GetPointer(0) + cellIdx * numPts;
    output->InsertNextCell(type, numPts, ptr);
  }

  return 1;
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::TraverseDualRecursively(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkHyperTreeGrid* input)
{
  // TODO - code duplication is evil, fold this with other TraverseDual
  // Create cell corner if cursor is at leaf
  if (cursor->IsLeaf())
  {
    // Center is a leaf, create dual items depending on dimension
    switch (input->GetDimension())
    {
      case 1:
        this->GenerateDualCornerFromLeaf1D(cursor, input);
        break;
      case 2:
        this->GenerateDualCornerFromLeaf2D(cursor, input);
        break;
      case 3:
        this->GenerateDualCornerFromLeaf3D(cursor, input);
        break;
    } // switch ( this->Dimension )
  }   // if ( cursor->IsLeaf() )
  else
  {
    // Cursor is not at leaf, recurse to all children
    int numChildren = input->GetNumberOfChildren();
    for (int child = 0; child < numChildren; ++child)
    {
      if (this->CheckAbort())
      {
        break;
      }
      cursor->ToChild(child);
      // Recurse
      this->TraverseDualRecursively(cursor, input);
      cursor->ToParent();
    } // child
  }   // else
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::GenerateDualCornerFromLeaf1D(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkHyperTreeGrid* input)
{
  // With d=1:
  //   (d-0)-faces are corners, neighbor cursors are 0 and 2
  //   (d-1)-faces do not exist
  //   (d-2)-faces do not exist

  // Retrieve neighbor (left/right) cursors
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorL =
    cursor->GetOrientedGeometryCursor(0);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorR =
    cursor->GetOrientedGeometryCursor(2);

  // Retrieve cursor center coordinates
  double pt[3];
  cursor->GetPoint(pt);

  // Check across d-face neighbors whether point must be adjusted
  if (!cursorL->HasTree())
  {
    // Move to left corner
    pt[input->GetOrientation()] -= .5 * cursor->GetSize()[input->GetOrientation()];
  }
  if (!cursorR->HasTree())
  {
    // Move to right corner
    pt[input->GetOrientation()] += .5 * cursor->GetSize()[input->GetOrientation()];
  }

  // Retrieve global index of center cursor
  vtkIdType id = cursor->GetGlobalNodeIndex();

  // Insert dual point at center of leaf cell
  this->Points->SetPoint(id, pt);

  // Storage for edge vertex IDs: dual cell ownership to cursor with higher index
  vtkIdType ids[2];
  ids[0] = id;

  // Check whether a dual edge to left neighbor exists
  if (cursorL->HasTree() && cursorL->IsLeaf())
  {
    // If left neighbor is a leaf, always create an edge
    ids[1] = cursorL->GetGlobalNodeIndex();
    this->Connectivity->InsertNextTypedTuple(ids);
  }

  // Check whether a dual edge to right neighbor exists
  if ((cursorR->HasTree() && cursorR->IsLeaf()) && cursorR->GetLevel() != cursor->GetLevel())
  {
    // If right neighbor is a leaf, create an edge only if right cell at higher level
    ids[1] = cursorR->GetGlobalNodeIndex();
    this->Connectivity->InsertNextTypedTuple(ids);
  }
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::GenerateDualCornerFromLeaf2D(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkHyperTreeGrid* input)
{
  // With d=2:
  //   (d-0)-faces are edges, neighbor cursors are 1, 3, 5, 7
  //   (d-1)-faces are corners, neighbor cursors are 0, 2, 6, 8
  //   (d-2)-faces do not exist

  // Retrieve (d-0)-neighbor (south/east/west/north) cursors
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorS =
    cursor->GetOrientedGeometryCursor(1);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorW =
    cursor->GetOrientedGeometryCursor(3);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorE =
    cursor->GetOrientedGeometryCursor(5);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorN =
    cursor->GetOrientedGeometryCursor(7);

  // Retrieve (d-1)-neighbor (southwest/southeast/northwest/northeast) cursors
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorSW =
    cursor->GetOrientedGeometryCursor(0);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorSE =
    cursor->GetOrientedGeometryCursor(2);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorNW =
    cursor->GetOrientedGeometryCursor(6);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorNE =
    cursor->GetOrientedGeometryCursor(8);

  // Retrieve 2D axes (east-west/south-north)
  unsigned int axisWE = input->GetOrientation() ? 0 : 1;
  unsigned int axisSN = input->GetOrientation() == 2 ? 1 : 2;

  // Retrieve cursor center coordinates
  double pt[3];
  cursor->GetPoint(pt);

  // Compute potential shifts
  double shift[2];
  shift[0] = .5 * cursor->GetSize()[axisWE];
  shift[1] = .5 * cursor->GetSize()[axisSN];

  // Check across edge neighbors whether point must be adjusted
  if (!cursorS->HasTree())
  {
    // Move to south edge
    pt[axisSN] -= shift[1];
  }
  if (!cursorW->HasTree())
  {
    // Move to west edge
    pt[axisWE] -= shift[0];
  }
  if (!cursorE->HasTree())
  {
    // Move to east edge
    pt[axisWE] += shift[0];
  }
  if (!cursorN->HasTree())
  {
    // Move to north edge
    pt[axisSN] += shift[1];
  }

  // Retrieve global index of center cursor
  vtkIdType id = cursor->GetGlobalNodeIndex();

  // Insert dual point at center of leaf cell
  this->Points->SetPoint(id, pt);

  // Storage for edge vertex IDs: dual cell ownership to cursor with higher index
  vtkIdType ids[4];
  ids[0] = id;

  // Retrieve current level to break corner ownership ties
  unsigned int level = cursor->GetLevel();

  // Check whether a dual cell around SW corner exists
  if (cursorSW->HasTree() && cursorSW->IsLeaf() && cursorS->HasTree() && cursorS->IsLeaf() &&
    cursorW->HasTree() && cursorW->IsLeaf())
  {
    // If SW, S, and W neighbors are leaves, always create a face
    ids[1] = cursorW->GetGlobalNodeIndex();
    ids[2] = cursorS->GetGlobalNodeIndex();
    ids[3] = cursorSW->GetGlobalNodeIndex();
    this->Connectivity->InsertNextTypedTuple(ids);
  }

  // Check whether a dual cell around SE corner exists
  if (cursorS->HasTree() && cursorS->IsLeaf() && cursorSE->HasTree() && cursorSE->IsLeaf() &&
    cursorE->HasTree() && cursorE->IsLeaf() && level != cursorE->GetLevel())
  {
    // If S, SE, and E neighbors are leaves, create a face if E at higher level
    ids[1] = cursorE->GetGlobalNodeIndex();
    ids[2] = cursorS->GetGlobalNodeIndex();
    ids[3] = cursorSE->GetGlobalNodeIndex();
    this->Connectivity->InsertNextTypedTuple(ids);
  }

  // Check whether a dual cell around NE corner exists
  if (cursorE->HasTree() && cursorE->IsLeaf() && cursorNE->HasTree() && cursorNE->IsLeaf() &&
    cursorN->HasTree() && cursorN->IsLeaf() && level != cursorE->GetLevel() &&
    level != cursorNE->GetLevel() && level != cursorN->GetLevel())
  {
    // If E, NE, and N neighbors are leaves, create a face if E, NE, N at higher level
    ids[1] = cursorE->GetGlobalNodeIndex();
    ids[2] = cursorN->GetGlobalNodeIndex();
    ids[3] = cursorNE->GetGlobalNodeIndex();
    this->Connectivity->InsertNextTypedTuple(ids);
  }

  // Check whether a dual cell around NW corner exists
  if (cursorW->HasTree() && cursorW->IsLeaf() && cursorN->HasTree() && cursorN->IsLeaf() &&
    cursorNW->HasTree() && cursorNW->IsLeaf() && level != cursorNW->GetLevel() &&
    level != cursorN->GetLevel())
  {
    // If W, N, and NW neighbors are leaves, create a face if NW and N at higher level
    ids[1] = cursorW->GetGlobalNodeIndex();
    ids[2] = cursorN->GetGlobalNodeIndex();
    ids[3] = cursorNW->GetGlobalNodeIndex();
    this->Connectivity->InsertNextTypedTuple(ids);
  }
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::GenerateDualCornerFromLeaf3D(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkHyperTreeGrid* vtkNotUsed(input))
{
  // With d=3:
  //   (d-0)-faces are faces, neighbor cursors are 4, 10, 12, 14, 16, 22
  //   (d-1)-faces are edges, neighbor cursors are 1, 3, 5, 7, 9, 11, 15, 17, 19, 21, 23, 25
  //   (d-2)-faces are corners, neighbor cursors are 0, 2, 6, 8, 18, 20, 24, 26

  // Retrieve cursors
  std::vector<vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor>> cursors;
  cursors.resize(27);
  for (unsigned int c = 0; c < 27; ++c)
  {
    cursors[c] = cursor->GetOrientedGeometryCursor(c);
  }

  // Retrieve cursor center coordinates
  double pt[3];
  cursor->GetPoint(pt);

  // Compute potential shifts
  double shift[3];
  shift[0] = .5 * cursor->GetSize()[0];
  shift[1] = .5 * cursor->GetSize()[1];
  shift[2] = .5 * cursor->GetSize()[2];

  // Index offset relative to center cursor (13)
  unsigned int offset = 1;

  // Check across face neighbors whether point must be adjusted
  for (unsigned int axis = 0; axis < 3; ++axis, offset *= 3)
  {
    if (!cursors[13 - offset]->HasTree())
    {
      // Move to negative side along axis
      pt[axis] -= shift[axis];
    }
    if (!cursors[13 + offset]->HasTree())
    {
      // Move to positive side along axis
      pt[axis] += shift[axis];
    }
  } // axis

  // Retrieve global index of center cursor
  vtkIdType id = cursor->GetGlobalNodeIndex();

  // Insert dual point at center of leaf cell
  this->Points->SetPoint(id, pt);

  // Storage for edge vertex IDs: dual cell ownership to cursor with higher index
  vtkIdType ids[8];

  // Retrieve current level to break corner ownership ties
  unsigned int level = cursor->GetLevel();

  // Iterate over leaf corners
  for (unsigned int c = 0; c < 8; ++c)
  {
    // Assume center cursor leaf owns the corner
    bool owner = true;

    // Iterate over every leaf touching the corner
    for (unsigned int l = 0; l < 8 && owner; ++l)
    {
      // Retrieve cursor index of touching leaf
      unsigned int index = CornerNeighborCursorsTable3D[c][l];

      // Compute whether corner is owned by another leaf
      if (index != 13)
      {
        if (!cursors[index]->HasTree() || !cursors[index]->IsLeaf() ||
          (cursors[index]->GetLevel() == level && index > 13))
        {
          // If neighbor leaf is out of bounds or has not been
          // refined to a leaf, this leaf does not own the corner
          // A level tie is broken in favor of the largest index
          owner = false;
        }
        else
        {
          // Collect the leaf indices for the dual cell
          ids[l] = cursors[index]->GetGlobalNodeIndex();
        }
      }
      else
      { // if ( index != 13 )
        // Collect the leaf indices for the dual cell
        ids[l] = cursors[index]->GetGlobalNodeIndex();
      } // else
    }   // l

    // If leaf owns the corner, create dual cell
    if (owner)
    {
      this->Connectivity->InsertNextTypedTuple(ids);
    } // if ( owner )
  }   // c
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::TraverseDualRecursively(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkBitArray* mask, vtkHyperTreeGrid* input)
{
  // TODO - code duplication is evil, fold this with other TraverseDual
  // Create cell corner if cursor is at leaf
  if (cursor->IsLeaf())
  {
    // Cursor is at leaf, retrieve its global index
    vtkIdType id = cursor->GetGlobalNodeIndex();

    // Center is a leaf, create dual items depending on dimension
    if (mask->GetValue(id))
    {
      switch (input->GetDimension())
      {
        case 2:
          this->ShiftDualCornerFromMaskedLeaf2D(cursor, mask, input);
          break;
        case 3:
          this->ShiftDualCornerFromMaskedLeaf3D(cursor, mask, input);
      } // switch ( this->Dimension )
    }
    else
    {
      switch (input->GetDimension())
      {
        case 1:
          this->GenerateDualCornerFromLeaf1D(cursor, input);
          break;
        case 2:
          this->GenerateDualCornerFromLeaf2D(cursor, mask, input);
          break;
        case 3:
          this->GenerateDualCornerFromLeaf3D(cursor, mask, input);
      } // switch ( this->Dimension )
    }   // else
  }     // if ( cursor->IsLeaf() )
  else
  {
    // Cursor is not at leaf, recurse to all children
    int numChildren = input->GetNumberOfChildren();
    for (int child = 0; child < numChildren; ++child)
    {
      if (this->CheckAbort())
      {
        break;
      }
      cursor->ToChild(child);
      // Recurse
      this->TraverseDualRecursively(cursor, mask, input);
      cursor->ToParent();
    } // child
  }   // else
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::ShiftDualCornerFromMaskedLeaf2D(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkBitArray* mask, vtkHyperTreeGrid* input)
{
  // With d=2:
  //   (d-0)-faces are edges, neighbor cursors are 1, 3, 5, 7
  //   (d-1)-faces are corners, neighbor cursors are 0, 2, 6, 8
  //   (d-2)-faces do not exist

  // Retrieve (d-0)-neighbor (south/east/west/north) cursors
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorS =
    cursor->GetOrientedGeometryCursor(1);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorW =
    cursor->GetOrientedGeometryCursor(3);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorE =
    cursor->GetOrientedGeometryCursor(5);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorN =
    cursor->GetOrientedGeometryCursor(7);

  // Retrieve (d-1)-neighbor (southwest/southeast/northwest/northeast) cursors
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorSW =
    cursor->GetOrientedGeometryCursor(0);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorSE =
    cursor->GetOrientedGeometryCursor(2);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorNW =
    cursor->GetOrientedGeometryCursor(6);
  vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorNE =
    cursor->GetOrientedGeometryCursor(8);

  // Retrieve global indices of non-center cursors

  // Retrieve 2D axes (east-west/south-north)
  unsigned int axisWE = input->GetOrientation() ? 0 : 1;
  unsigned int axisSN = input->GetOrientation() == 2 ? 1 : 2;

  // Retrieve current level to break corner ownership ties
  unsigned int level = cursor->GetLevel();

  // Check whether dual point across S edge must be adjusted
  if (cursorS->HasTree() && cursorS->IsLeaf() && cursorS->GetLevel() < level)
  {
    vtkIdType idS = cursorS->GetGlobalNodeIndex();
    if (!mask->GetValue(idS))
    {
      // Dual point must be adjusted
      this->PointShifted[idS] = true;
      this->PointShifts[axisSN][idS] =
        cursorS->GetTree()->GetScale(axisSN) * this->ReductionFactors[cursorS->GetLevel()];
    }
  }

  // Check whether dual point across W edge must be adjusted
  if (cursorW->HasTree() && cursorW->IsLeaf() && cursorW->GetLevel() < level)
  {
    vtkIdType idW = cursorW->GetGlobalNodeIndex();
    if (!mask->GetValue(idW))
    {
      // Dual point must be adjusted
      this->PointShifted[idW] = true;
      this->PointShifts[axisWE][idW] =
        cursorW->GetTree()->GetScale(axisWE) * this->ReductionFactors[cursorW->GetLevel()];
    }
  }

  // Check whether dual point across E face must be adjusted
  if (cursorE->HasTree() && cursorE->IsLeaf() && cursorE->GetLevel() < level)
  {
    vtkIdType idE = cursorE->GetGlobalNodeIndex();
    if (!mask->GetValue(idE))
    {
      // Dual point must be adjusted
      this->PointShifted[idE] = true;
      this->PointShifts[axisWE][idE] =
        -cursorE->GetTree()->GetScale(axisWE) * this->ReductionFactors[cursorE->GetLevel()];
    }
  }

  // Check whether dual point across N edge must be adjusted
  if (cursorN->HasTree() && cursorN->IsLeaf() && cursorN->GetLevel() < level)
  {
    vtkIdType idN = cursorN->GetGlobalNodeIndex();
    if (!mask->GetValue(idN))
    {
      // Dual point must be adjusted
      this->PointShifted[idN] = true;
      this->PointShifts[axisSN][idN] =
        -cursorN->GetTree()->GetScale(axisSN) * this->ReductionFactors[cursorN->GetLevel()];
    }
  }

  // Check whether dual point across SE corner must be adjusted
  if (cursorSE->HasTree() && cursorSE->IsLeaf() && cursorSE->GetLevel() < level)
  {
    vtkIdType idSE = cursorSE->GetGlobalNodeIndex();
    if (!mask->GetValue(idSE) && !this->PointShifted[idSE])
    {
      // Dual point must be adjusted
      double shift[3];
      cursorSE->GetTree()->GetScale(shift);
      double factor = this->ReductionFactors[cursorSE->GetLevel()];
      this->PointShifts[axisWE][idSE] = factor * shift[axisWE];
      this->PointShifts[axisSN][idSE] = factor * shift[axisSN];
    }
  }

  // Check whether dual point across SW corner must be adjusted
  if (cursorSW->HasTree() && cursorSW->IsLeaf() && cursorSW->GetLevel() < level)
  {
    vtkIdType idSW = cursorSW->GetGlobalNodeIndex();
    if (!mask->GetValue(idSW) && !this->PointShifted[idSW])
    {
      // Dual point must be adjusted
      double shift[3];
      cursorSW->GetTree()->GetScale(shift);
      double factor = this->ReductionFactors[cursorSW->GetLevel()];
      this->PointShifts[axisWE][idSW] = -factor * shift[axisWE];
      this->PointShifts[axisSN][idSW] = factor * shift[axisSN];
    }
  }

  // Check whether dual point across NW corner must be adjusted
  if (cursorNW->HasTree() && cursorNW->IsLeaf() && cursorNW->GetLevel() < level)
  {
    vtkIdType idNW = cursorNW->GetGlobalNodeIndex();
    if (!mask->GetValue(idNW) && !this->PointShifted[idNW])
    {
      // Dual point must be adjusted
      double shift[3];
      cursorNW->GetTree()->GetScale(shift);
      double factor = this->ReductionFactors[cursorNW->GetLevel()];
      this->PointShifts[axisWE][idNW] = factor * shift[axisWE];
      this->PointShifts[axisSN][idNW] = -factor * shift[axisSN];
    }
  }

  // Check whether dual point across NE corner must be adjusted
  if (cursorNE->HasTree() && cursorNE->IsLeaf() && cursorNE->GetLevel() < level)
  {
    vtkIdType idNE = cursorNE->GetGlobalNodeIndex();
    if (!mask->GetValue(idNE) && !this->PointShifted[idNE])
    {
      // Dual point must be adjusted
      double shift[3];
      cursorNE->GetTree()->GetScale(shift);
      double factor = this->ReductionFactors[cursorNE->GetLevel()];
      this->PointShifts[axisWE][idNE] = -factor * shift[axisWE];
      this->PointShifts[axisSN][idNE] = -factor * shift[axisSN];
    }
  }
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::ShiftDualCornerFromMaskedLeaf3D(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkBitArray* mask,
  vtkHyperTreeGrid* vtkNotUsed(input))
{
  // With d=3:
  //   (d-0)-faces are faces, neighbor cursors are 4, 10, 12, 14, 16, 22
  //   (d-1)-faces are edges, neighbor cursors are 1, 3, 5, 7, 9, 11, 15, 17, 19, 21, 23, 25
  //   (d-2)-faces are corners, neighbor cursors are 0, 2, 6, 8, 18, 20, 24, 26

  // Retrieve current level to break corner ownership ties
  unsigned int level = cursor->GetLevel();

  // Check whether dual points across face neighbors must be adjusted
  int offset = 1;
  for (unsigned int axis = 0; axis < 3; ++axis, offset *= 3)
  {
    vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorM =
      cursor->GetOrientedGeometryCursor(13 - offset);
    if (cursorM->HasTree() && cursorM->IsLeaf() && cursorM->GetLevel() < level)
    {
      vtkIdType idM = cursorM->GetGlobalNodeIndex();
      if (!mask->GetValue(idM))
      {
        // Dual point must be adjusted
        this->PointShifted[idM] = true;
        this->PointShifts[axis][idM] =
          cursorM->GetTree()->GetScale(axis) * this->ReductionFactors[cursorM->GetLevel()];
      }
    }
    vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorP =
      cursor->GetOrientedGeometryCursor(13 + offset);
    if (cursorP->HasTree() && cursorP->IsLeaf() && cursorP->GetLevel() < level)
    {
      vtkIdType idP = cursorP->GetGlobalNodeIndex();
      if (!mask->GetValue(idP))
      {
        // Dual point must be adjusted
        this->PointShifted[idP] = true;
        this->PointShifts[axis][idP] =
          -cursorP->GetTree()->GetScale(axis) * this->ReductionFactors[cursorP->GetLevel()];
      }
    }
  } // axis

  // Check whether dual points across edge neighbors must be adjusted
  int i = 1;
  for (int axis1 = 0; axis1 < 2; ++axis1, i *= 3)
  {
    int j = 3 * i;
    for (int axis2 = axis1 + 1; axis2 < 3; ++axis2, j *= 3)
    {
      for (int o2 = -1; o2 < 2; o2 += 2)
      {
        for (int o1 = -1; o1 < 2; o1 += 2)
        {
          int index = 13 + o1 * (i * o2 + j);
          vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorE =
            cursor->GetOrientedGeometryCursor(static_cast<unsigned int>(index));
          if (cursorE->HasTree() && cursorE->IsLeaf() && cursorE->GetLevel() < level)
          {
            vtkIdType idE = cursorE->GetGlobalNodeIndex();
            if (!mask->GetValue(idE) && !this->PointShifted[idE])
            {
              // Dual point must be adjusted
              this->PointShifted[idE] = true;
              double shift[3];
              cursorE->GetTree()->GetScale(shift);
              double factor = this->ReductionFactors[cursorE->GetLevel()];
              this->PointShifts[axis1][idE] = -o1 * o2 * factor * shift[axis1];
              this->PointShifts[axis2][idE] = -o1 * factor * shift[axis2];
            }
          }
        } // o1
      }   // o2
    }     // axis2
  }       // axis1

  // Check whether dual points across corner neighbors must be adjusted
  for (int o3 = -1; o3 < 2; o3 += 2)
  {
    for (int o2 = -1; o2 < 2; o2 += 2)
    {
      offset = o2 * (o3 + 3) + 9;
      for (int o1 = -1; o1 < 2; o1 += 2)
      {
        int index = 13 + o1 * offset;
        vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorC =
          cursor->GetOrientedGeometryCursor(index);
        if (cursorC->HasTree() && cursorC->IsLeaf() && cursorC->GetLevel() < level)
        {
          vtkIdType idC = cursorC->GetGlobalNodeIndex();
          if (!mask->GetValue(idC) && !this->PointShifted[idC])
          {
            // Dual point must be adjusted
            this->PointShifted[idC] = true;
            double shift[3];
            cursorC->GetTree()->GetScale(shift);
            double factor = this->ReductionFactors[cursorC->GetLevel()];
            this->PointShifts[0][idC] = -o1 * o2 * o3 * factor * shift[0];
            this->PointShifts[1][idC] = -o1 * o2 * factor * shift[1];
            this->PointShifts[2][idC] = -o1 * factor * shift[2];
          }
        }
      } // o1
    }   // o2
  }     // o3
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::GenerateDualCornerFromLeaf2D(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkBitArray* mask, vtkHyperTreeGrid* input)
{
  // With d=2:
  //   (d-0)-faces are edges, neighbor cursors are 1, 3, 5, 7
  //   (d-1)-faces are corners, neighbor cursors are 0, 2, 6, 8
  //   (d-2)-faces do not exist

  const unsigned int dS = 1;
  const unsigned int dW = 3;
  const unsigned int dE = 5;
  const unsigned int dN = 7;
  const unsigned int dSW = 0;
  const unsigned int dSE = 2;
  const unsigned int dNW = 6;
  const unsigned int dNE = 8;

  // Retrieve 2D axes (east-west/south-north)
  unsigned int axisWE = input->GetOrientation() ? 0 : 1;
  unsigned int axisSN = input->GetOrientation() == 2 ? 1 : 2;

  // Retrieve cursor center coordinates
  double pt[3];
  cursor->GetPoint(pt);

  // Compute potential shifts
  double shift[2];
  shift[0] = .5 * cursor->GetSize()[axisWE];
  shift[1] = .5 * cursor->GetSize()[axisSN];

  // When a mask is present, edge as well as face shifts are possible
  bool shifted = false;

  // Check across edge neighbors whether point must be adjusted
  if (!cursor->HasTree(dS))
  {
    // Move to south edge
    pt[axisSN] -= shift[1];
    shifted = true;
  }
  else if (cursor->IsLeaf(dS) && mask->GetValue(cursor->GetGlobalNodeIndex(dS)))
  {
    // Move to south edge
    pt[axisSN] -= shift[1];
    shifted = true;
  }

  if (!cursor->HasTree(dW))
  {
    // Move to west edge
    pt[axisWE] -= shift[0];
    shifted = true;
  }
  else if (cursor->IsLeaf(dW) && mask->GetValue(cursor->GetGlobalNodeIndex(dW)))
  {
    // Move to west edge
    pt[axisWE] -= shift[0];
    shifted = true;
  }

  if (!cursor->HasTree(dE))
  {
    // Move to east edge
    pt[axisWE] += shift[0];
    shifted = true;
  }
  else if (cursor->IsLeaf(dE) && mask->GetValue(cursor->GetGlobalNodeIndex(dE)))
  {
    // Move to east edge
    pt[axisWE] += shift[0];
    shifted = true;
  }

  if (!cursor->HasTree(dN))
  {
    // Move to north edge
    pt[axisSN] += shift[1];
    shifted = true;
  }
  else if (cursor->IsLeaf(dN) && mask->GetValue(cursor->GetGlobalNodeIndex(dN)))
  {
    // Move to north edge
    pt[axisSN] += shift[1];
    shifted = true;
  }

  // Only when point was not moved to edge, check corner neighbors
  if (!shifted)
  {
    if (!cursor->HasTree(dSW))
    {
      // Move to southwest corner
      pt[axisWE] -= shift[0];
      pt[axisSN] -= shift[1];
    }
    else if (cursor->IsLeaf(dSW) && mask->GetValue(cursor->GetGlobalNodeIndex(dSW)))
    {
      // Move to southwest corner
      pt[axisWE] -= shift[0];
      pt[axisSN] -= shift[1];
    }

    if (!cursor->HasTree(dSE))
    {
      // Move to southeast corner
      pt[axisWE] += shift[0];
      pt[axisSN] -= shift[1];
    }
    else if (cursor->IsLeaf(dSE) && mask->GetValue(cursor->GetGlobalNodeIndex(dSE)))
    {
      // Move to southeast corner
      pt[axisWE] += shift[0];
      pt[axisSN] -= shift[1];
    }

    if (!cursor->HasTree(dNW))
    {
      // Move to northwest corner
      pt[axisWE] -= shift[0];
      pt[axisSN] += shift[1];
    }
    else if (cursor->IsLeaf(dNW) && mask->GetValue(cursor->GetGlobalNodeIndex(dNW)))
    {
      // Move to northwest corner
      pt[axisWE] -= shift[0];
      pt[axisSN] += shift[1];
    }

    if (!cursor->HasTree(dNE))
    {
      // Move to northeast corner
      pt[axisWE] += shift[0];
      pt[axisSN] += shift[1];
    }
    else if (cursor->IsLeaf(dNE) && mask->GetValue(cursor->GetGlobalNodeIndex(dNE)))
    {
      // Move to northeast corner
      pt[axisWE] += shift[0];
      pt[axisSN] += shift[1];
    }
  } // if ( ! shifted )

  // Retrieve global index of center cursor
  vtkIdType id = cursor->GetGlobalNodeIndex();

  // Insert dual point at center of leaf cell
  assert(id < input->GetNumberOfCells());
  this->Points->SetPoint(id, pt);

  // If cell is masked, terminate recursion, no dual cell will be generated
  if (mask->GetValue(id))
  {
    return;
  }

  // Storage for edge vertex IDs: dual cell ownership to cursor with higher index
  vtkIdType ids[4];
  ids[0] = id;

  // Retrieve current level to break corner ownership ties
  unsigned int level = cursor->GetLevel();

  // Check whether a dual cell around SW corner exists
  if (cursor->HasTree(dSW) && cursor->HasTree(dS) && cursor->HasTree(dW) && cursor->IsLeaf(dSW) &&
    cursor->IsLeaf(dS) && cursor->IsLeaf(dW))
  {
    vtkIdType idSW, idS, idW;
    if (!mask->GetValue(idSW = cursor->GetGlobalNodeIndex(dSW)) &&
      !mask->GetValue(idS = cursor->GetGlobalNodeIndex(dS)) &&
      !mask->GetValue(idW = cursor->GetGlobalNodeIndex(dW)))
    {
      // If SW, S, and W neighbors are leaves, always create a face
      ids[1] = idW;
      ids[2] = idS;
      ids[3] = idSW;
      this->Connectivity->InsertNextTypedTuple(ids);
    }
  }

  // Check whether a dual cell around SE corner exists
  if (cursor->HasTree(dS) && cursor->HasTree(dSE) && cursor->HasTree(dE) && cursor->IsLeaf(dS) &&
    cursor->IsLeaf(dSE) && cursor->IsLeaf(dE))
  {
    vtkIdType idS, idSE, idE;
    if (!mask->GetValue(idS = cursor->GetGlobalNodeIndex(dS)) &&
      !mask->GetValue(idSE = cursor->GetGlobalNodeIndex(dSE)) &&
      !mask->GetValue(idE = cursor->GetGlobalNodeIndex(dE)) && level != cursor->GetLevel(dE))
    {
      // If S, SE, and E neighbors are leaves, create a face if E at higher level
      ids[1] = idE;
      ids[2] = idS;
      ids[3] = idSE;
      this->Connectivity->InsertNextTypedTuple(ids);
    }
  }

  // Check whether a dual cell around NE corner exists
  if (cursor->HasTree(dE) && cursor->HasTree(dNE) && cursor->HasTree(dN) && cursor->IsLeaf(dE) &&
    cursor->IsLeaf(dNE) && cursor->IsLeaf(dN))
  {
    vtkIdType idE, idNE, idN;
    if (!mask->GetValue(idE = cursor->GetGlobalNodeIndex(dE)) &&
      !mask->GetValue(idNE = cursor->GetGlobalNodeIndex(dNE)) &&
      !mask->GetValue(idN = cursor->GetGlobalNodeIndex(dN)) && level != cursor->GetLevel(dE) &&
      level != cursor->GetLevel(dNE) && level != cursor->GetLevel(dN))
    {
      // If E, NE, and N neighbors are leaves, create a face if E, NE, N at higher level
      ids[1] = idE;
      ids[2] = idN;
      ids[3] = idNE;
      this->Connectivity->InsertNextTypedTuple(ids);
    }
  }

  // Check whether a dual cell around NW corner exists
  if (cursor->HasTree(dW) && cursor->HasTree(dN) && cursor->HasTree(dNW) && cursor->IsLeaf(dW) &&
    cursor->IsLeaf(dN) && cursor->IsLeaf(dNW))
  {
    vtkIdType idW, idN, idNW;
    if (!mask->GetValue(idW = cursor->GetGlobalNodeIndex(dW)) &&
      !mask->GetValue(idN = cursor->GetGlobalNodeIndex(dN)) &&
      !mask->GetValue(idNW = cursor->GetGlobalNodeIndex(dNW)) && level != cursor->GetLevel(dNW) &&
      level != cursor->GetLevel(dN))
    {
      // If W, N, and NW neighbors are leaves, create a face if NW and N at higher level
      ids[1] = idW;
      ids[2] = idN;
      ids[3] = idNW;
      this->Connectivity->InsertNextTypedTuple(ids);
    }
  }
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridToDualGrid::GenerateDualCornerFromLeaf3D(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, vtkBitArray* mask,
  vtkHyperTreeGrid* vtkNotUsed(input))
{
  // With d=3:
  //   (d-0)-faces are faces, neighbor cursors are 4, 10, 12, 14, 16, 22
  //   (d-1)-faces are edges, neighbor cursors are 1, 3, 5, 7, 9, 11, 15, 17, 19, 21, 23, 25
  //   (d-2)-faces are corners, neighbor cursors are 0, 2, 6, 8, 18, 20, 24, 26

  // Retrieve cursor center coordinates
  double pt[3];
  cursor->GetPoint(pt);

  // Compute potential shifts
  double shift[3];
  shift[0] = .5 * cursor->GetSize()[0];
  shift[1] = .5 * cursor->GetSize()[1];
  shift[2] = .5 * cursor->GetSize()[2];

  // When a mask is present, corner, edge, and face shifts are possible
  bool shifted = false;

  // Check across face neighbors whether point must be adjusted
  unsigned int offset = 1;
  for (unsigned int axis = 0; axis < 3; ++axis, offset *= 3)
  {
    vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorM =
      cursor->GetOrientedGeometryCursor(13 - offset);
    if (!cursorM->HasTree())
    {
      // Move to negative side along axis
      pt[axis] -= shift[axis];
      shifted = true;
    }
    else
    {
      vtkIdType idM = cursorM->GetGlobalNodeIndex();
      if (cursorM->IsLeaf() && mask->GetValue(idM))
      {
        // Move to negative side along axis
        pt[axis] -= shift[axis];
        shifted = true;
      }
    }
    vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorP =
      cursor->GetOrientedGeometryCursor(13 + offset);
    if (!cursorP->HasTree())
    {
      // Move to positive side along axis
      pt[axis] += shift[axis];
      shifted = true;
    }
    else
    {
      vtkIdType idP = cursorP->GetGlobalNodeIndex();
      if (cursorP->IsLeaf() && mask->GetValue(idP))
      {
        // Move to positive side along axis
        pt[axis] += shift[axis];
        shifted = true;
      }
    }
  } // axis

  // Only when point was not moved to face, check edge neighbors
  if (!shifted)
  {
    int i = 1;
    for (int axis1 = 0; axis1 < 2; ++axis1, i *= 3)
    {
      int j = 3 * i;
      for (int axis2 = axis1 + 1; axis2 < 3; ++axis2, j *= 3)
      {
        for (int o2 = -1; o2 < 2; o2 += 2)
        {
          for (int o1 = -1; o1 < 2; o1 += 2)
          {
            int index = 13 + o1 * (i * o2 + j);
            vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorE =
              cursor->GetOrientedGeometryCursor(static_cast<unsigned int>(index));
            if (!cursorE->HasTree())
            {
              // Move to corresponding edge
              pt[axis1] += o1 * o2 * shift[axis1];
              pt[axis2] += o1 * shift[axis2];
              shifted = true;
            }
            else
            {
              vtkIdType idE = cursorE->GetGlobalNodeIndex();
              if (cursorE->IsLeaf() && mask->GetValue(idE))
              {
                // Move to corresponding edge
                pt[axis1] += o1 * o2 * shift[axis1];
                pt[axis2] += o1 * shift[axis2];
                shifted = true;
              }
            }
          } // o1
        }   // o2
      }     // axis2
    }       // axis1
  }         // if ( ! shifted )

  // Only when point was neither moved to face nor to edge, check corners neighbors
  if (!shifted)
  {
    // Iterate over all 8 corners
    for (int o3 = -1; o3 < 2; o3 += 2)
    {
      for (int o2 = -1; o2 < 2; o2 += 2)
      {
        offset = o2 * (o3 + 3) + 9;
        for (int o1 = -1; o1 < 2; o1 += 2)
        {
          int index = 13 + o1 * static_cast<int>(offset);
          vtkSmartPointer<vtkHyperTreeGridOrientedGeometryCursor> cursorC =
            cursor->GetOrientedGeometryCursor(static_cast<unsigned int>(index));
          if (!cursorC->HasTree())
          {
            // Move to corresponding corner
            pt[0] += o1 * o2 * o3 * shift[0];
            pt[1] += o1 * o2 * shift[1];
            pt[2] += o1 * shift[2];
          }
          else
          { // if cursor
            vtkIdType idC = cursorC->GetGlobalNodeIndex();
            if (cursorC->IsLeaf() && mask->GetValue(idC))
            {
              // Move to corresponding corner
              pt[0] += o1 * o2 * o3 * shift[0];
              pt[1] += o1 * o2 * shift[1];
              pt[2] += o1 * shift[2];
            }
          } // if cursor
        }   // o1
      }     // o2
    }       // o3
  }         // if ( ! shifted )

  // Retrieve global index of center cursor
  vtkIdType id = cursor->GetGlobalNodeIndex();

  // Insert dual point at center of leaf cell
  // assert ( id < input->GetGlobalNodeIndexMax() + 1 );
  this->Points->SetPoint(id, pt);

  // Storage for edge vertex IDs: dual cell ownership to cursor with higher index
  vtkIdType ids[8];

  // Retrieve current level to break corner ownership ties
  unsigned int level = cursor->GetLevel();

  // Iterate over leaf corners
  for (unsigned int c = 0; c < 8; ++c)
  {
    // Assume center cursor leaf owns the corner
    bool owner = true;
    unsigned int real_l = 0;

    // Iterate over every leaf touching the corner
    for (unsigned int l = 0; l < 8 && owner; ++l)
    {
      // Retrieve cursor index of touching leaf
      unsigned int index = CornerNeighborCursorsTable3D[c][l];

      // Compute whether corner is owned by another leaf
      if (index != 13)
      {
        if (!cursor->HasTree(index) || !cursor->IsLeaf(index) ||
          (cursor->GetLevel(index) == level && index > 13))
        {
          // If neighbor leaf is out of bounds or has not been
          // refined to a leaf, this leaf does not own the corner
          // A level tie is broken in favor of the largest index
          owner = false;
        }
        else
        {
          vtkIdType idglobal = cursor->GetGlobalNodeIndex(index);
          if (!mask->GetValue(idglobal))
          {
            // Collect the leaf indices for the dual cell
            ids[real_l++] = cursor->GetGlobalNodeIndex(index);
          }
          else
          {
            owner = false;
          }
        }
      }
      else
      { // if ( index != 13 )
        // Collect the leaf indices for the dual cell
        ids[real_l++] = id;
      } // else
    }   // l

    // If leaf owns the corner, create dual cell
    if (owner)
    {
      if (real_l != 8)
      {
        if (real_l == 0)
        {
          continue;
        }
        vtkIdType last = ids[real_l - 1];
        for (; real_l < 8; ++real_l)
        {
          ids[real_l] = last;
        }
      }
      this->Connectivity->InsertNextTypedTuple(ids);
    } // if ( owner )
  }   // c
}
VTK_ABI_NAMESPACE_END