File: vtkHyperTreeGridPlaneCutter.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 (883 lines) | stat: -rw-r--r-- 25,245 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
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
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#include "vtkHyperTreeGridPlaneCutter.h"

#include "vtkBitArray.h"
#include "vtkCellArray.h"
#include "vtkCellData.h"
#include "vtkCleanPolyData.h"
#include "vtkCutter.h"
#include "vtkDataObject.h"
#include "vtkHyperTreeGrid.h"
#include "vtkHyperTreeGridNonOrientedGeometryCursor.h"
#include "vtkHyperTreeGridNonOrientedMooreSuperCursor.h"
#include "vtkIdList.h"
#include "vtkInformation.h"
#include "vtkMath.h"
#include "vtkObjectFactory.h"
#include "vtkPlane.h"
#include "vtkPointData.h"
#include "vtkPoints.h"
#include "vtkPolyData.h"
#include "vtkUnstructuredGrid.h"

#include <cassert>
#include <cmath>

VTK_ABI_NAMESPACE_BEGIN
namespace
{
vtkIdType First8Integers[] = { 0, 1, 2, 3, 4, 5, 6, 7 };

constexpr unsigned int MooreCursors3D[26] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 14, 15, 16,
  17, 18, 19, 20, 21, 22, 23, 24, 25, 26 };
}

vtkStandardNewMacro(vtkHyperTreeGridPlaneCutter);

//------------------------------------------------------------------------------
vtkHyperTreeGridPlaneCutter::vtkHyperTreeGridPlaneCutter()
{
  this->Points = nullptr;
  this->Cells = nullptr;

  // Initialize plane parameters
  std::fill(this->Plane, this->Plane + 4, 0.);

  // By default a non-conforming output mesh is produced for better rendering
  this->Dual = 0;

  // By default member variables for dual-based computation are not used
  this->SelectedCells = nullptr;
  this->Centers = nullptr;
  this->Cutter = nullptr;
  this->Leaves = nullptr;
}

//------------------------------------------------------------------------------
vtkHyperTreeGridPlaneCutter::~vtkHyperTreeGridPlaneCutter()
{
  if (this->Points)
  {
    this->Points->Delete();
    this->Points = nullptr;
  }

  if (this->Cells)
  {
    this->Cells->Delete();
    this->Cells = nullptr;
  }

  if (this->Leaves)
  {
    this->Leaves->Delete();
    this->Leaves = nullptr;
  }

  if (this->Centers)
  {
    this->Centers->Delete();
    this->Centers = nullptr;
  }

  if (this->Cutter)
  {
    this->Cutter->Delete();
    this->Cutter = nullptr;
  }

  if (this->SelectedCells)
  {
    this->SelectedCells->Delete();
    this->SelectedCells = nullptr;
  }
}

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

  os << indent << "Plane: ( " << this->Plane[0] << " ) * X + ( " << this->Plane[1] << " ) * Y + ( "
     << this->Plane[2] << " ) * Z = " << this->Plane[3] << "\n";

  if (this->Dual)
  {
    os << indent << "Dual: Yes\n";
  }
  else
  {
    os << indent << "Dual: No\n";
  }

  if (this->Points)
  {
    os << indent << "Points:\n";
    this->Points->PrintSelf(os, indent.GetNextIndent());
  }
  else
  {
    os << indent << "Points: ( none )\n";
  }

  if (this->Cells)
  {
    os << indent << "Cells:\n";
    this->Cells->PrintSelf(os, indent.GetNextIndent());
  }
  else
  {
    os << indent << "Cells: ( none )\n";
  }

  if (this->Leaves)
  {
    os << indent << "Leaves:\n";
    this->Leaves->PrintSelf(os, indent.GetNextIndent());
  }
  else
  {
    os << indent << "Leaves: ( none )\n";
  }

  if (this->Centers)
  {
    os << indent << "Centers:\n";
    this->Centers->PrintSelf(os, indent.GetNextIndent());
  }
  else
  {
    os << indent << "Centers: ( none )\n";
  }

  if (this->Cutter)
  {
    os << indent << "Cutter:\n";
    this->Cutter->PrintSelf(os, indent.GetNextIndent());
  }
  else
  {
    os << indent << "Cutter: ( none )\n";
  }
}

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

//------------------------------------------------------------------------------
void vtkHyperTreeGridPlaneCutter::Reset()
{
  // Points and Cells are created in the constructor
  if (this->Points)
  {
    this->Points->Delete();
  }
  this->Points = vtkPoints::New();
  if (this->Cells)
  {
    this->Cells->Delete();
  }
  this->Cells = vtkCellArray::New();
  if (this->Centers)
  {
    this->Centers->Initialize();
  }
  if (this->Leaves)
  {
    this->Leaves->Initialize();
  }
  if (this->Cutter)
  {
    this->Cutter->SetNumberOfContours(0);
  }
  if (this->SelectedCells)
  {
    this->SelectedCells->Initialize();
  }
}

//------------------------------------------------------------------------------
int vtkHyperTreeGridPlaneCutter::ProcessTrees(vtkHyperTreeGrid* input, vtkDataObject* outputDO)
{
  vtkPolyData* output = vtkPolyData::SafeDownCast(outputDO);
  if (!output)
  {
    vtkErrorMacro("Incorrect type of output: " << outputDO->GetClassName());
    return 0;
  }

  // This filter works only with 3D grids
  if (input->GetDimension() != 3)
  {
    vtkErrorMacro(<< "Bad input dimension:" << input->GetDimension());
    return 0;
  }

  // Reset Data
  this->Reset();
  output->Initialize();

  // Retrieve input point data
  this->InData = input->GetCellData();

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

  // Compute cut on dual or primal input depending on specification
  if (this->Dual)
  {
    // Initialize output point data
    this->OutData = output->GetPointData();
    this->OutData->CopyAllocate(this->InData);

    // Storage for leaf indices
    if (this->Leaves == nullptr)
    {
      this->Leaves = vtkIdList::New();
    }
    this->Leaves->SetNumberOfIds(8);

    // Initialize storage for dual geometry
    if (this->Centers == nullptr)
    {
      this->Centers = vtkPoints::New();
    }
    this->Centers->SetNumberOfPoints(8);

    // Convert plane parameters into normal/origin specification
    unsigned int maxId = 0;
    if (fabs(this->Plane[1]) > fabs(this->Plane[0]))
    {
      maxId = 1;
    }
    if (fabs(this->Plane[2]) > fabs(this->Plane[maxId]))
    {
      maxId = 2;
    }
    double origin[] = { 0., 0., 0. };
    origin[maxId] = this->Plane[3] / this->Plane[maxId];
    vtkPlane* plane = vtkPlane::New();
    plane->SetOrigin(origin);
    plane->SetNormal(this->Plane[0], this->Plane[1], this->Plane[2]);

    // Initialize plane cutter
    if (this->Cutter == nullptr)
    {
      this->Cutter = vtkCutter::New();
    }
    this->Cutter->GenerateTrianglesOff();
    this->Cutter->SetCutFunction(plane);

    // Clean up
    plane->Delete();

    // Create storage to keep track of selected cells
    if (this->SelectedCells == nullptr)
    {
      this->SelectedCells = vtkBitArray::New();
    }
    vtkIdType numCells = input->GetNumberOfCells();
    this->SelectedCells->SetNumberOfTuples(numCells);
    for (vtkIdType i = 0; i < numCells; ++i)
    {
      // Initialization is needed because not all cells are pre-processed
      this->SelectedCells->SetValue(i, 0);
    }

    // First pass across tree roots to evince cells intersected by contours
    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 input tree
      input->InitializeNonOrientedGeometryCursor(cursor, index);
      // Pre-process tree recursively
      this->RecursivelyPreProcessTree(cursor);
    } // it

    // Second pass across tree roots: now compute isocontours recursively
    input->InitializeTreeIterator(it);
    vtkNew<vtkHyperTreeGridNonOrientedMooreSuperCursor> supercursor;
    while (it.GetNextTree(index))
    {
      if (this->CheckAbort())
      {
        break;
      }
      // Initialize new Moore cursor at root of current tree
      input->InitializeNonOrientedMooreSuperCursor(supercursor, index);
      // Generate leaf cell centers recursively
      this->RecursivelyProcessTreeDual(supercursor);
    } // it

    // Clean up
    this->SelectedCells->Delete();
    this->SelectedCells = nullptr;
  } // if ( this->Dual )
  else
  {
    // Initialize output cell data
    this->OutData = output->GetCellData();
    this->OutData->CopyAllocate(this->InData);

    // 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);
      // Generate leaf cell centers recursively
      this->RecursivelyProcessTreePrimal(cursor);
    } // it
  }   // else

  // Set output geometry and topology
  output->SetPoints(this->Points);
  this->Points->FastDelete();
  this->Points = nullptr;
  output->SetPolys(this->Cells);
  this->Cells->FastDelete();
  this->Cells = nullptr;

  // Clean and squeeze output
  vtkCleanPolyData* cleaner = vtkCleanPolyData::New();
  cleaner->ConvertPolysToLinesOff();
  cleaner->SetInputData(output);
  cleaner->Update();
  output->ShallowCopy(cleaner->GetOutput());
  output->Squeeze();
  cleaner->Delete();
  return 1;
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridPlaneCutter::RecursivelyProcessTreePrimal(
  vtkHyperTreeGridNonOrientedGeometryCursor* cursor)
{
  // If cursor is at a masked cell stop recursion
  vtkIdType inId = cursor->GetGlobalNodeIndex();
  if (this->InMask && this->InMask->GetValue(inId))
  {
    return;
  }

  // Retrieve cursor geometry
  double* origin = cursor->GetOrigin();
  double* size = cursor->GetSize();

  // Initialize cell coordinates
  double cellCoords[8][3];
  for (int i = 0; i < 8; ++i)
  {
    cellCoords[i][0] = (i & 1) ? origin[0] + size[0] : origin[0];
    // Checking if the plane is equal to the boundary of a cell.
    // If it is, we need to shift it a tiny bit.
    // Check is done on all axis.
    // NOTE: we set cellCoords to std::sqrt(VTK_DBL_MIN) if the plane passes by the origin
    // because distance computation, needed later, requires squaring those values.
    // Since VTK_DBL_MIN is the smallest normal double value, VTK_DBL_MIN*VTK_DBL_MIN == 0,
    // and sqrt(VTK_DBL_MIN)*std::sqrt(VTK_DBL_MIN) == VTK_DBL_MIN, which is what we want
    if (this->IsPlaneOrthogonalToXAxis())
    {
      if (cellCoords[i][0] == this->Plane[3])
      {
        cellCoords[i][0] += std::abs(cellCoords[i][0]) > std::sqrt(VTK_DBL_MIN)
          ? VTK_DBL_EPSILON * std::abs(cellCoords[i][0])
          : std::sqrt(VTK_DBL_MIN);
      }
    }
    cellCoords[i][1] = (i & 2) ? origin[1] + size[1] : origin[1];
    if (this->IsPlaneOrthogonalToYAxis())
    {
      if (cellCoords[i][1] == this->Plane[3])
      {
        cellCoords[i][1] += std::abs(cellCoords[i][1]) > std::sqrt(VTK_DBL_MIN)
          ? VTK_DBL_EPSILON * std::abs(cellCoords[i][1])
          : std::sqrt(VTK_DBL_MIN);
      }
    }
    cellCoords[i][2] = (i & 4) ? origin[2] + size[2] : origin[2];
    if (this->IsPlaneOrthogonalToZAxis())
    {
      if (cellCoords[i][2] == this->Plane[3])
      {
        cellCoords[i][2] += std::abs(cellCoords[i][2]) > std::sqrt(VTK_DBL_MIN)
          ? VTK_DBL_EPSILON * std::abs(cellCoords[i][2])
          : std::sqrt(VTK_DBL_MIN);
      }
    }
  }

  const double SQRT_DBL_EPSILON = std::sqrt(VTK_DBL_EPSILON);

  // Check cell-plane intersection
  double functEval[8];
  if (this->CheckIntersection(cellCoords, functEval))
  {
    // Create plane cut if cursor is at leaf
    if (cursor->IsLeaf())
    {
      // Keep track of the number of intersection points
      int n = 0;

      // Storage for intersection points
      double points[8][3];

      // Iterate over cell vertices
      for (int i = 0; i < 8; ++i)
      {
        // Check all cell edges
        if (std::abs(functEval[i]) < SQRT_DBL_EPSILON)
        {
          // If current vertex is intersected then save it
          memcpy(points[n], cellCoords[i], 3 * sizeof(double));
          ++n;
        }
        else
        {
          // Check every edge of the current vertex.
          if (!(i & 1) && functEval[i] * functEval[i + 1] < 0)
          {
            // Edge in X
            this->PlaneCut(i, i + 1, cellCoords, n, points);
          }
          if (!(i & 2) && functEval[i] * functEval[i + 2] < 0)
          {
            // Edge in Y
            this->PlaneCut(i, i + 2, cellCoords, n, points);
          }
          if (!(i & 4) && functEval[i] * functEval[i + 4] < 0)
          {
            // Edge in Z
            this->PlaneCut(i, i + 4, cellCoords, n, points);
          }
        } // else
      }   // i

      // Now reorder points if necessary
      this->ReorderCutPoints(n, points);

      // Storage for face vertex IDs
      vtkIdType ids[8];
      for (int i = 0; i < n; ++i)
      {
        // Save points and get their IDs
        ids[i] = this->Points->InsertNextPoint(points[i]);
      }

      // Insert next face
      vtkIdType outId = this->Cells->InsertNextCell(n, ids);

      // Copy face data from that of the cell from which it comes
      this->OutData->CopyData(this->InData, inId, outId);
    } // 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->RecursivelyProcessTreePrimal(cursor);
        cursor->ToParent();
      } // ichild
    }   // else
  }     // CheckIntersection
}

//------------------------------------------------------------------------------
bool vtkHyperTreeGridPlaneCutter::RecursivelyPreProcessTree(
  vtkHyperTreeGridNonOrientedGeometryCursor* cursor)
{
  // If cursor is at a masked cell stop recursion
  vtkIdType id = cursor->GetGlobalNodeIndex();
  if (this->InMask && this->InMask->GetValue(id))
  {
    return false;
  }

  // A node is not selected until proven otherwise
  bool selected = false;

  // Retrieve cursor geometry
  double* origin = cursor->GetOrigin();
  double* size = cursor->GetSize();

  // Initialize cell coordinates
  double cellCoords[8][3];
  for (int i = 0; i < 8; ++i)
  {
    cellCoords[i][0] = (i & 1) ? origin[0] + size[0] : origin[0];
    cellCoords[i][1] = (i & 2) ? origin[1] + size[1] : origin[1];
    cellCoords[i][2] = (i & 4) ? origin[2] + size[2] : origin[2];
  }

  // Check cell-plane intersection
  if (this->CheckIntersection(cellCoords))
  {
    // Selected this node
    if (cursor->IsLeaf())
    {
      selected = true;
    } // 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 and keep track of whether this branch is selected
        selected |= this->RecursivelyPreProcessTree(cursor);
        cursor->ToParent();
      } // ichild
    }   // else
  }     // if ( this->CheckIntersection )

  // Update list of selected cells
  this->SelectedCells->SetTuple1(id, selected);

  // Return whether current node was selected
  return selected;
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridPlaneCutter::RecursivelyProcessTreeDual(
  vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor)
{
  // If cursor is at a masked cell stop recursion
  vtkIdType id = cursor->GetGlobalNodeIndex();
  if (this->InMask && this->InMask->GetValue(id))
  {
    return;
  }

  // Descend further into input trees only if cursor is not a leaf
  if (!cursor->IsLeaf())
  {
    // Check if cursor is at selected cell
    if (!this->SelectedCells->GetTuple1(id))
    {
      // Cell is not selected until proven otherwise
      bool selected = false;

      // Iterate over all cursors of Von Neumann neighborhood around center
      for (unsigned int neighbor = 0; neighbor < 26 && !selected; ++neighbor)
      {
        // Retrieve global index of neighbor
        unsigned int indN = MooreCursors3D[neighbor];
        if (cursor->HasTree(indN))
        {
          vtkIdType idN = cursor->GetGlobalNodeIndex(indN);

          // Decide whether neighbor was selected
          selected = (this->SelectedCells->GetTuple1(idN) != 0.0);
        }
        else
        {
          selected = false;
        }
      } // neighbor

      // No dual cell with a corner at cursor center will be intersected
      if (!selected)
      {
        return;
      }
    } // if ( this->SelectedCells->GetTuple1( id ) )

    // Recurse to all children
    int numChildren = cursor->GetNumberOfChildren();
    for (int ichild = 0; ichild < numChildren; ++ichild)
    {
      if (this->CheckAbort())
      {
        break;
      }
      cursor->ToChild(ichild);
      // Recurse
      this->RecursivelyProcessTreeDual(cursor);
      cursor->ToParent();
    } // ichild
  }   // if ( ! cursor->IsLeaf() )
  else
  {
    // Cursor is at leaf, iterate over its corners
    for (unsigned int cornerIdx = 0; cornerIdx < 8; ++cornerIdx)
    {
      if (this->CheckAbort())
      {
        break;
      }
      // Cell is not selected until proven otherwise
      bool owner = true;

      // Iterate over every leaf touching the corner and check ownership
      for (unsigned int leafIdx = 0; leafIdx < 8 && owner; ++leafIdx)
      {
        owner = cursor->GetCornerCursors(cornerIdx, leafIdx, this->Leaves);
      } // leafIdx

      // If cell owns dual cell, compute intersection thereof
      if (owner)
      {
        // Create dual cell to be intersected as unstructured grid
        vtkUnstructuredGrid* dual = vtkUnstructuredGrid::New();
        dual->Allocate(1, 1);
        dual->InsertNextCell(VTK_VOXEL, 8, First8Integers);
        dual->GetPointData()->CopyAllocate(this->InData);

        // Iterate over cell corners
        double x[3];
        for (int _cornerIdx = 0; _cornerIdx < 8; ++_cornerIdx)
        {
          // Get cursor corresponding to this corner
          vtkIdType cursorId = this->Leaves->GetId(_cornerIdx);

          // Retrieve neighbor coordinates and store them
          cursor->GetPoint(cursorId, x);
          this->Centers->SetPoint(_cornerIdx, x);

          // Retrieve neighbor index and corresponding input scalar value
          vtkIdType idN = cursor->GetGlobalNodeIndex(cursorId);

          // Assign scalar value attached to this cell corner
          dual->GetPointData()->CopyData(this->InData, idN, _cornerIdx);
        } // _cornerIdx

        // Assign geometry of dual cell
        dual->SetPoints(this->Centers);

        // Compute intersection with plane
        this->Cutter->SetInputData(dual);
        this->Cutter->Update();

        // Append computed polygons if some are present in cutter output
        vtkPolyData* pd = this->Cutter->GetOutput();
        vtkIdType nPoints = pd->GetNumberOfPoints();
        if (nPoints)
        {
          // Keep handle to cut point data
          vtkPointData* pdata = pd->GetPointData();

          // Append new points to existing cut points
          vtkIdType offset = this->Points->GetNumberOfPoints();
          double pt[3];
          for (vtkIdType i = 0; i < nPoints; ++i)
          {
            // Retrieve cut point coordinates and insert them into output points
            pd->GetPoint(i, pt);
            this->Points->InsertNextPoint(pt);

            // Copy cut point data to that of corresponding output point
            this->OutData->CopyData(pdata, i, i + offset);
          } // i

          // Append new elements to existing cut element
          vtkIdType ids[8];
          vtkIdType nCells = pd->GetNumberOfCells();
          for (int i = 0; i < nCells; ++i)
          {
            // Retrieve element vertex ids
            vtkIdList* vertices = pd->GetCell(i)->GetPointIds();

            // Appropriately offset vertex ids
            vtkIdType n = vertices->GetNumberOfIds();
            for (int j = 0; j < n; ++j)
            {
              ids[j] = vertices->GetId(j) + offset;
            } // j

            // Insert next cell with offset ids
            this->Cells->InsertNextCell(n, ids);
          } // i
        }   // if ( nPoints )

        // Clean up
        dual->Delete();
      } // if ( owner )
    }   // cornerIdx
  }     // else
}

//------------------------------------------------------------------------------
bool vtkHyperTreeGridPlaneCutter::CheckIntersection(double cellCoords[8][3], double functEval[8])
{
  // Iterate over cell vertices
  int i;
  for (i = 0; i < 8; ++i)
  {
    // Evaluate the plane in every coordinate
    functEval[i] = cellCoords[i][0] * this->Plane[0] + cellCoords[i][1] * this->Plane[1] +
      cellCoords[i][2] * this->Plane[2] - this->Plane[3];
  } // i

  // Evaluate plane equation at first corner
  double firstVal = functEval[0];

  // Check if there is any sign change
  i = 7;
  while (i && functEval[i] * firstVal > 0.)
  {
    --i;
  }

  // Intersection if while statement broke early
  return (i != 0);
}

//------------------------------------------------------------------------------
bool vtkHyperTreeGridPlaneCutter::CheckIntersection(double cellCoords[8][3])
{
  // Evaluate plane equation at first corner
  double firstVal = cellCoords[0][0] * this->Plane[0] + cellCoords[0][1] * this->Plane[1] +
    cellCoords[0][2] * this->Plane[2] - this->Plane[3];

  // Check if there is any sign change
  int i = 1;
  bool sameSign = true;
  while (i < 8 && sameSign)
  {
    double functEval = cellCoords[i][0] * this->Plane[0] + cellCoords[i][1] * this->Plane[1] +
      cellCoords[i][2] * this->Plane[2] - this->Plane[3];
    sameSign = (firstVal * functEval > 0);
    ++i;
  }

  return !sameSign;
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridPlaneCutter::SetPlane(double a, double b, double c, double d)
{
  assert(!(a == 0 && b == 0 && c == 0) && "Plane's normal equals zero");
  this->Plane[0] = a;
  this->Plane[1] = b;
  this->Plane[2] = c;
  this->Plane[3] = d;
  if (a == 0.0 && b == 0.0)
  {
    this->AxisAlignment = 2;
  }
  else if (b == 0.0 && c == 0.0)
  {
    this->AxisAlignment = 0;
  }
  else if (a == 0.0 && c == 0.0)
  {
    this->AxisAlignment = 1;
  }
  else
  {
    this->AxisAlignment = -1;
  }
  this->Modified();
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridPlaneCutter::PlaneCut(
  int i, int j, double cellCoords[8][3], int& n, double point[][3])
{
  // The distance between vertex indices gives us the axis direction
  if (j - i == 1)
  {
    // X direction
    point[n][0] =
      (this->Plane[3] - this->Plane[1] * cellCoords[i][1] - this->Plane[2] * cellCoords[i][2]) /
      this->Plane[0];
    point[n][1] = cellCoords[i][1];
    point[n][2] = cellCoords[i][2];
  } // if ( j - i == 1 )
  else if (j - i == 2)
  {
    // Y direction
    point[n][0] = cellCoords[i][0];
    point[n][1] =
      (this->Plane[3] - this->Plane[0] * cellCoords[i][0] - this->Plane[2] * cellCoords[i][2]) /
      this->Plane[1];
    point[n][2] = cellCoords[i][2];
  }    // else if ( j - i == 2 )
  else // if ( j - i == 4 )
  {
    // Z direction
    point[n][0] = cellCoords[i][0];
    point[n][1] = cellCoords[i][1];
    point[n][2] =
      (this->Plane[3] - this->Plane[0] * cellCoords[i][0] - this->Plane[1] * cellCoords[i][1]) /
      this->Plane[2];
  } // else if ( j - i == 4 )

  // Move to next point
  ++n;
}

//------------------------------------------------------------------------------
void vtkHyperTreeGridPlaneCutter::ReorderCutPoints(int n, double points[][3])
{
  // Iterate over all polygonal vertices but the last one
  for (int i = 0; i < n - 2; ++i)
  {
    // Search the closest point to i in the sense of sharing the most coordinates
    int index = i + 1; // in practice a don't care, set to satisfy -Wuninitialized
    int minDistance = 4;
    for (int j = i + 1; j < n; ++j)
    {
      // Compute the distance: number of different coordinates
      int distance = 0;
      if (points[j][0] != points[i][0])
      {
        ++distance;
      }
      if (points[j][1] != points[i][1])
      {
        ++distance;
      }
      if (points[j][2] != points[i][2])
      {
        ++distance;
      }
      if (distance < minDistance)
      {
        // Store the index of current point and update minDistance
        index = j;
        minDistance = distance;
      }
    }
    if (index != i + 1)
    {
      // If the closest point is not i + 1, then swap point positions
      double swap[3];
      memcpy(swap, points[index], 3 * sizeof(double));
      memcpy(points[index], points[i + 1], 3 * sizeof(double));
      memcpy(points[i + 1], swap, 3 * sizeof(double));
    }
  }
}
VTK_ABI_NAMESPACE_END