vtkVolumeOfRevolutionFilter.cxx 18.6 KB
Newer Older
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
/*=========================================================================

  Program:   Visualization Toolkit
  Module:    vtkVolumeOfRevolutionFilter.cxx

  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
  All rights reserved.
  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.

     This software is distributed WITHOUT ANY WARRANTY; without even
     the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
     PURPOSE.  See the above copyright notice for more information.

=========================================================================*/
#include "vtkVolumeOfRevolutionFilter.h"

#include "vtkCellArray.h"
#include "vtkCellData.h"
#include "vtkCellIterator.h"
#include "vtkGenericCell.h"
#include "vtkIdList.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkMath.h"
#include "vtkMergePoints.h"
#include "vtkNew.h"
#include "vtkObjectFactory.h"
#include "vtkPointData.h"
#include "vtkPolyData.h"
#include "vtkUnstructuredGrid.h"

#include <vector>

vtkStandardNewMacro(vtkVolumeOfRevolutionFilter);

namespace
{
struct AxisOfRevolution
{
  double Position[3];
  double Direction[3];
};

void RevolvePoint(const double in[3], const AxisOfRevolution* axis,
                  double angleInRadians, double out[3])
{
  double c = cos(angleInRadians);
  double cm = 1. - c;
  double s  = sin(angleInRadians);

  double translated[3];
  vtkMath::Subtract(in, axis->Position, translated);

  double dot = vtkMath::Dot(translated, axis->Direction);
  double cross[3];
  vtkMath::Cross(translated, axis->Direction, cross);

  for (vtkIdType i=0;i<3;i++)
59
  {
60
61
    out[i] = ((translated[i]*c + axis->Direction[i]*dot*cm - cross[i]*s) +
              axis->Position[i]);
62
  }
63
64
65
66
67
68
69
70
71
72
73
74
75
}

void RevolvePoints(vtkDataSet* pts, vtkPoints* newPts, AxisOfRevolution* axis,
                   double sweepAngle, int resolution, vtkPointData *outPd,
                   bool partialSweep)
{
  double angleInRadians = vtkMath::RadiansFromDegrees(sweepAngle/resolution);

  vtkIdType n2DPoints = pts->GetNumberOfPoints();
  vtkIdType counter = 0;
  double p2d[3], p3d[3];

  for (int i=0; i<resolution + partialSweep; i++)
76
  {
77
    for (int id = 0; id < n2DPoints; id++)
78
    {
79
80
81
82
83
84
      pts->GetPoint(id, p2d);
      RevolvePoint(p2d, axis, i * angleInRadians, p3d);
      newPts->SetPoint(counter, p3d);
      outPd->CopyData(pts->GetPointData(), i, counter);
      counter++;
    }
85
  }
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
}

template <int CellType>
void Revolve(vtkIdList* pointIds, vtkIdType n2DPoints, int resolution,
             vtkCellArray *connectivity, vtkUnsignedCharArray *types,
             vtkIdTypeArray *locations, vtkCellData *inCd, vtkIdType cellId,
             vtkCellData *outCd, bool partialSweep);

template <>
void Revolve<VTK_VERTEX>(vtkIdList* pointIds, vtkIdType n2DPoints,
                         int resolution, vtkCellArray *connectivity,
                         vtkUnsignedCharArray *types, vtkIdTypeArray *locations,
                         vtkCellData *inCd, vtkIdType cellId,
                         vtkCellData *outCd, bool partialSweep)
{
  vtkIdType newPtIds[2], newCellId;

  newPtIds[0] = pointIds->GetId(0);

  for (int i=0; i<resolution; i++)
106
  {
107
108
109
    newPtIds[1] = (pointIds->GetId(0) +
                   ((i+1)%(resolution + partialSweep))*n2DPoints);
    newCellId = connectivity->InsertNextCell(2, newPtIds);
110
    locations->InsertNextValue(connectivity->GetInsertLocation(0));
111
112
113
    types->InsertNextValue(VTK_LINE);
    outCd->CopyData(inCd, cellId, newCellId);
    newPtIds[0] = newPtIds[1];
114
  }
115
116
117
118
119
120
121
122
123
124
125
126
127
}

template <>
void Revolve<VTK_POLY_VERTEX>(vtkIdList* pointIds, vtkIdType n2DPoints,
                              int resolution, vtkCellArray *connectivity,
                              vtkUnsignedCharArray *types,
                              vtkIdTypeArray *locations, vtkCellData *inCd,
                              vtkIdType cellId, vtkCellData *outCd,
                              bool partialSweep)
{
  vtkNew<vtkIdList> pointId;
  pointId->SetNumberOfIds(1);
  for (vtkIdType i=0; i<pointIds->GetNumberOfIds(); i++)
128
  {
129
130
131
132
    pointId->SetId(0,pointIds->GetId(i));
    Revolve<VTK_VERTEX>(pointId.GetPointer(), n2DPoints, resolution,
                        connectivity, types, locations, inCd, cellId, outCd,
                        partialSweep);
133
  }
134
135
136
137
138
139
140
141
142
143
144
145
146
}

template <>
void Revolve<VTK_LINE>(vtkIdList* pointIds, vtkIdType n2DPoints, int resolution,
                       vtkCellArray *connectivity, vtkUnsignedCharArray *types,
                       vtkIdTypeArray *locations, vtkCellData *inCd,
                       vtkIdType cellId, vtkCellData *outCd, bool partialSweep)
{
  static const int nPoints = 2;

  vtkIdType newPtIds[2*nPoints], newCellId;

  for (vtkIdType i=0;i<nPoints;i++)
147
  {
148
    newPtIds[i] = pointIds->GetId(i);
149
  }
150
151

  for (int i=0; i<resolution; i++)
152
  {
153
    for (vtkIdType j=0; j<nPoints; j++)
154
    {
155
156
      newPtIds[2*nPoints-1-j] =
        pointIds->GetId(j) + ((i+1)%(resolution + partialSweep))*n2DPoints;
157
    }
158
    newCellId = connectivity->InsertNextCell(2*nPoints, newPtIds);
159
    locations->InsertNextValue(connectivity->GetInsertLocation(0));
160
161
162
    types->InsertNextValue(VTK_QUAD);
    outCd->CopyData(inCd, cellId, newCellId);
    for (vtkIdType j=0;j<nPoints;j++)
163
    {
164
165
      newPtIds[nPoints-1-j] = newPtIds[j+nPoints];
    }
166
  }
167
168
169
170
171
172
173
174
175
176
177
178
179
180
}

template <>
void Revolve<VTK_POLY_LINE>(vtkIdList* pointIds, vtkIdType n2DPoints,
                            int resolution, vtkCellArray *connectivity,
                            vtkUnsignedCharArray *types,
                            vtkIdTypeArray *locations, vtkCellData *inCd,
                            vtkIdType cellId, vtkCellData *outCd,
                            bool partialSweep)
{
  vtkNew<vtkIdList> newPointIds;
  newPointIds->SetNumberOfIds(2);
  newPointIds->SetId(0,pointIds->GetId(0));
  for (vtkIdType i=1; i<pointIds->GetNumberOfIds(); i++)
181
  {
182
183
184
185
186
    newPointIds->SetId(1,pointIds->GetId(i));
    Revolve<VTK_LINE>(newPointIds.GetPointer(), n2DPoints, resolution,
                      connectivity, types, locations, inCd, cellId, outCd,
                      partialSweep);
    newPointIds->SetId(0,pointIds->GetId(i));
187
  }
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
}

template <>
void Revolve<VTK_TRIANGLE>(vtkIdList* pointIds, vtkIdType n2DPoints,
                           int resolution, vtkCellArray *connectivity,
                           vtkUnsignedCharArray *types,
                           vtkIdTypeArray *locations, vtkCellData *inCd,
                           vtkIdType cellId, vtkCellData *outCd,
                           bool partialSweep)
{
  static const int nPoints = 3;

  vtkIdType newPtIds[2*nPoints], newCellId;

  for (vtkIdType i=0;i<nPoints;i++)
203
  {
204
    newPtIds[i] = pointIds->GetId(i);
205
  }
206
207

  for (int i=0; i<resolution; i++)
208
  {
209
    for (vtkIdType j=0; j<nPoints; j++)
210
    {
211
212
      newPtIds[j+nPoints] =
        pointIds->GetId(j) + ((i+1)%(resolution + partialSweep))*n2DPoints;
213
    }
214
    newCellId = connectivity->InsertNextCell(2*nPoints, newPtIds);
215
    locations->InsertNextValue(connectivity->GetInsertLocation(0));
216
217
218
    types->InsertNextValue(VTK_WEDGE);
    outCd->CopyData(inCd, cellId, newCellId);
    for (vtkIdType j=0;j<nPoints;j++)
219
    {
220
221
      newPtIds[j] = newPtIds[j+nPoints];
    }
222
  }
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
}

template <>
void Revolve<VTK_TRIANGLE_STRIP>(vtkIdList* pointIds, vtkIdType n2DPoints,
                                 int resolution, vtkCellArray *connectivity,
                                 vtkUnsignedCharArray *types,
                                 vtkIdTypeArray *locations, vtkCellData *inCd,
                                 vtkIdType cellId, vtkCellData *outCd,
                                 bool partialSweep)
{
  vtkNew<vtkIdList> newPointIds;
  newPointIds->SetNumberOfIds(3);
  newPointIds->SetId(0,pointIds->GetId(0));
  newPointIds->SetId(1,pointIds->GetId(1));
  for (vtkIdType i=2; i<pointIds->GetNumberOfIds(); i++)
238
  {
239
240
241
242
243
244
    newPointIds->SetId(2,pointIds->GetId(i));
    Revolve<VTK_TRIANGLE>(newPointIds.GetPointer(), n2DPoints, resolution,
                      connectivity, types, locations, inCd, cellId, outCd,
                          partialSweep);
    newPointIds->SetId(0,pointIds->GetId(i));
    newPointIds->SetId(1,pointIds->GetId(i-1));
245
  }
246
247
248
249
250
251
252
253
254
255
256
257
258
}

template <>
void Revolve<VTK_QUAD>(vtkIdList* pointIds, vtkIdType n2DPoints, int resolution,
                       vtkCellArray *connectivity, vtkUnsignedCharArray *types,
                       vtkIdTypeArray *locations, vtkCellData *inCd,
                       vtkIdType cellId, vtkCellData *outCd, bool partialSweep)
{
  static const int nPoints = 4;

  vtkIdType newPtIds[2*nPoints], newCellId;

  for (vtkIdType i=0;i<nPoints;i++)
259
  {
260
    newPtIds[i] = pointIds->GetId(i);
261
  }
262
263

  for (int i=0; i<resolution; i++)
264
  {
265
    for (vtkIdType j=0; j<nPoints; j++)
266
    {
267
268
      newPtIds[j+nPoints] =
        pointIds->GetId(j) + ((i+1)%(resolution + partialSweep))*n2DPoints;
269
    }
270
    newCellId = connectivity->InsertNextCell(2*nPoints, newPtIds);
271
    locations->InsertNextValue(connectivity->GetInsertLocation(0));
272
273
274
    types->InsertNextValue(VTK_HEXAHEDRON);
    outCd->CopyData(inCd, cellId, newCellId);
    for (vtkIdType j=0; j<nPoints; j++)
275
    {
276
277
      newPtIds[j] = newPtIds[j+nPoints];
    }
278
  }
279
280
281
282
283
284
285
286
287
288
289
290
291
292
}

template <>
void Revolve<VTK_PIXEL>(vtkIdList* pointIds, vtkIdType n2DPoints,
                        int resolution, vtkCellArray *connectivity,
                        vtkUnsignedCharArray *types, vtkIdTypeArray *locations,
                        vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
                        bool partialSweep)
{
  static const int nPoints = 4;

  vtkIdType newPtIds[2*nPoints], newCellId;

  for (vtkIdType i=0;i<nPoints;i++)
293
  {
294
    newPtIds[i] = pointIds->GetId(i);
295
  }
296
297

  for (int i=0; i<resolution; i++)
298
  {
299
    for (vtkIdType j=0; j<nPoints; j++)
300
    {
301
302
      newPtIds[j+nPoints] =
        pointIds->GetId(j) + ((i+1)%(resolution + partialSweep))*n2DPoints;
303
    }
304
    newCellId = connectivity->InsertNextCell(2*nPoints, newPtIds);
305
    locations->InsertNextValue(connectivity->GetInsertLocation(0));
306
307
308
    types->InsertNextValue(VTK_HEXAHEDRON);
    outCd->CopyData(inCd, cellId, newCellId);
    for (vtkIdType j=0; j<nPoints; j++)
309
    {
310
311
      newPtIds[j] = newPtIds[j+nPoints];
    }
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
}

template <>
void Revolve<VTK_POLYGON>(vtkIdList* pointIds, vtkIdType n2DPoints,
                          int resolution, vtkCellArray *connectivity,
                          vtkUnsignedCharArray *types,
                          vtkIdTypeArray *locations, vtkCellData *inCd,
                          vtkIdType cellId, vtkCellData *outCd,
                          bool partialSweep)
{
  // A swept polygon creates a polyhedron with two polygon faces and <nPoly>
  // quad faces, comprised from 2*<nPoly> points. Because polyhedra have a
  // special connectivity format, the length of the connectivity array is
  // 1 + (<nPoly>+2) + 2*<nPoly> + 4*<nPoly> = 7*<nPoly> + 3.
  // ^        ^           ^           ^
  // integer describing # of faces (<nPoly> + 2)
  //          ^           ^           ^
  //          integers describing # of points per face
  //                      ^           ^
  //                      point ids for the two polygon faces
  //                                  ^
  //                                  point ids for the 4 quad faces

  int nPoly = pointIds->GetNumberOfIds();
  std::vector<vtkIdType> newPtIds(7*nPoly + 3);
  // newFacePtIds are pointers to the point arrays describing each face
  std::vector<vtkIdType*> newFacePtIds(nPoly + 2);
  vtkIdType newCellId;

  newPtIds[0] = nPoly + 2;
  newPtIds[1] = nPoly;
  newPtIds[nPoly + 2] = nPoly;
  newFacePtIds[0] = &newPtIds[2];
  newFacePtIds[1] = &newPtIds[nPoly+3];
  for (vtkIdType i=0;i<nPoly;i++)
348
  {
349
350
351
352
    // All of the subsequent faces have four point ids
    newPtIds[3 + 2*nPoly + 5*i] = 4;
    newFacePtIds[2+i] = &newPtIds[4 + 2*nPoly + 5*i];
    newFacePtIds[0][i] = pointIds->GetId(i);
353
  }
354
355

  for (int i=0; i<resolution; i++)
356
  {
357
    for (vtkIdType j=0; j<nPoly; j++)
358
    {
359
360
      newFacePtIds[1][nPoly-1-j] =
        pointIds->GetId(j) + ((i+1)%(resolution + partialSweep))*n2DPoints;
361
    }
362
    for (vtkIdType j=0; j<nPoly; j++)
363
    {
364
365
366
367
      newFacePtIds[j+2][0] = newFacePtIds[0][j];
      newFacePtIds[j+2][1] = newFacePtIds[0][(j+1)%nPoly];
      newFacePtIds[j+2][2] = newFacePtIds[1][(2*nPoly-2-j)%nPoly];
      newFacePtIds[j+2][3] = newFacePtIds[1][nPoly-1-j];
368
    }
369
    newCellId = connectivity->InsertNextCell(7*nPoly+3, &newPtIds[0]);
370
    locations->InsertNextValue(connectivity->GetInsertLocation(0));
371
372
373
    types->InsertNextValue(VTK_POLYHEDRON);
    outCd->CopyData(inCd, cellId, newCellId);
    for (vtkIdType j=0; j<nPoly; j++)
374
    {
375
376
      newFacePtIds[0][j] = newFacePtIds[1][nPoly-1-j];
    }
377
  }
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
}

int RevolveCell(int cellType, vtkIdList* pointIds, vtkIdType n2DPoints,
                int resolution, vtkCellArray *connectivity,
                vtkUnsignedCharArray *types, vtkIdTypeArray *locations,
                vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
                bool partialSweep)
{
  int returnValue = 0;
#define RevolveCellCase(CellType)                                       \
  case CellType:                                                        \
    Revolve<CellType>(pointIds, n2DPoints, resolution, connectivity,    \
                      types, locations, inCd, cellId, outCd, partialSweep); \
    break

  switch (cellType)
394
  {
395
396
397
398
399
400
401
402
403
404
405
    RevolveCellCase(VTK_VERTEX);
    RevolveCellCase(VTK_POLY_VERTEX);
    RevolveCellCase(VTK_LINE);
    RevolveCellCase(VTK_POLY_LINE);
    RevolveCellCase(VTK_TRIANGLE);
    RevolveCellCase(VTK_TRIANGLE_STRIP);
    RevolveCellCase(VTK_POLYGON);
    RevolveCellCase(VTK_PIXEL);
    RevolveCellCase(VTK_QUAD);
    default:
      returnValue = 1;
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
  return returnValue;
#undef RevolveCellCase
}
}

//----------------------------------------------------------------------------
vtkVolumeOfRevolutionFilter::vtkVolumeOfRevolutionFilter()
{
  this->SweepAngle = 360.0;
  this->Resolution = 12; // 30 degree increments
  this->AxisPosition[0] = this->AxisPosition[1] = this->AxisPosition[2] = 0.;
  this->AxisDirection[0] = this->AxisDirection[1] = 0.;
  this->AxisDirection[2] = 1.;
  this->OutputPointsPrecision = DEFAULT_PRECISION;
}

//----------------------------------------------------------------------------
vtkVolumeOfRevolutionFilter::~vtkVolumeOfRevolutionFilter()
{
}

//----------------------------------------------------------------------------
int vtkVolumeOfRevolutionFilter::RequestData(
  vtkInformation *vtkNotUsed(request),
  vtkInformationVector **inputVector,
  vtkInformationVector *outputVector)
{
  // get the info objects
  vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
  vtkInformation *outInfo = outputVector->GetInformationObject(0);

  // get the input and output
  vtkDataSet *input = vtkDataSet::SafeDownCast(
    inInfo->Get(vtkDataObject::DATA_OBJECT()));
  vtkUnstructuredGrid *output = vtkUnstructuredGrid::SafeDownCast(
    outInfo->Get(vtkDataObject::DATA_OBJECT()));
  vtkPointData* inPd = input->GetPointData();
  vtkCellData* inCd = input->GetCellData();
  vtkPointData* outPd = output->GetPointData();
  vtkCellData* outCd = output->GetCellData();

  vtkNew<vtkPoints> outPts;

  // Check to see that the input data is amenable to this operation
  vtkCellIterator* it = input->NewCellIterator();
  for (it->InitTraversal(); !it->IsDoneWithTraversal(); it->GoToNextCell())
453
  {
454
455
    int cellDimension = it->GetCellDimension();
    if (cellDimension > 2)
456
    {
457
458
459
      vtkErrorMacro(<<"All cells must have a topological dimension < 2.");
      return 1;
    }
460
  }
461
462
463
464
  it->Delete();

  // Set up output points
  if (this->OutputPointsPrecision == vtkAlgorithm::DEFAULT_PRECISION)
465
  {
466
467
    vtkPointSet *inputPointSet = vtkPointSet::SafeDownCast(input);
    if (inputPointSet)
468
    {
469
      outPts->SetDataType(inputPointSet->GetPoints()->GetDataType());
470
    }
471
    else
472
    {
473
474
      outPts->SetDataType(VTK_FLOAT);
    }
475
  }
476
  else if (this->OutputPointsPrecision == vtkAlgorithm::SINGLE_PRECISION)
477
  {
478
    outPts->SetDataType(VTK_FLOAT);
479
  }
480
  else if (this->OutputPointsPrecision == vtkAlgorithm::DOUBLE_PRECISION)
481
  {
482
    outPts->SetDataType(VTK_DOUBLE);
483
  }
484
485
486
487

  // determine whether or not the sweep angle is a full 2*pi
  bool partialSweep = false;
  if (fabs(360. - fabs(this->SweepAngle)) > 1024*VTK_DBL_EPSILON)
488
  {
489
     partialSweep = true;
490
  }
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507

  // Set up output points and point data
  outPts->SetNumberOfPoints(input->GetNumberOfPoints() * (this->Resolution +
                                                          partialSweep));
  outPd->CopyAllocate(inPd, input->GetNumberOfPoints() * (this->Resolution +
                                                          partialSweep));

  // Set up output cell data
  vtkIdType nNewCells = input->GetNumberOfCells() * this->Resolution;
  outCd->CopyAllocate(inCd, nNewCells);

  vtkNew<vtkUnsignedCharArray> outTypes;
  vtkNew<vtkIdTypeArray> outLocations;
  vtkNew<vtkCellArray> outCells;

  AxisOfRevolution axis;
  for (vtkIdType i=0; i<3; i++)
508
  {
509
510
    axis.Position[i] = this->AxisPosition[i];
    axis.Direction[i] = this->AxisDirection[i];
511
  }
512
513
514
515
516
517

  RevolvePoints(input, outPts.GetPointer(), &axis, this->SweepAngle,
                this->Resolution, outPd, partialSweep);

  it = input->NewCellIterator();
  for (it->InitTraversal(); !it->IsDoneWithTraversal(); it->GoToNextCell())
518
  {
519
520
521
522
523
    if (RevolveCell(it->GetCellType(), it->GetPointIds(),
                    input->GetNumberOfPoints(), this->Resolution,
                    outCells.GetPointer(), outTypes.GetPointer(),
                    outLocations.GetPointer(), inCd, it->GetCellId(), outCd,
                    partialSweep) == 1)
524
    {
525
526
527
      vtkWarningMacro(<<"No method for revolving cell type "
                      << it->GetCellType() <<". Skipping.");
    }
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
  it->Delete();

  output->SetPoints(outPts.GetPointer());
  output->SetCells(outTypes.GetPointer(), outLocations.GetPointer(),
                   outCells.GetPointer());

  return 1;
}

//----------------------------------------------------------------------------
int vtkVolumeOfRevolutionFilter::FillInputPortInformation(
  int, vtkInformation *info)
{
  info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
  return 1;
}

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

  os << indent << "Resolution: " << this->Resolution << "\n";
  os << indent << "Sweep Angle: " << this->SweepAngle << "\n";
  os << indent << "Axis Position: (" << this->AxisPosition[0] << ","
     << this->AxisPosition[1] << "," << this->AxisPosition[2] << ")\n";
  os << indent << "Axis Direction: (" << this->AxisPosition[0] << ","
     << this->AxisDirection[1] << "," << this->AxisDirection[2] << ")\n";
  os << indent << "Output Points Precision: "<<this->OutputPointsPrecision
     << "\n";
}