TTK
Loading...
Searching...
No Matches
ttkIcosphere.cpp
Go to the documentation of this file.
1#include <ttkIcosphere.h>
2
3#include <ttkUtils.h>
4
5#include <vtkInformation.h>
6#include <vtkObjectFactory.h>
7
8#include <vtkCellArray.h>
9#include <vtkDoubleArray.h>
10#include <vtkFloatArray.h>
11#include <vtkIdTypeArray.h>
12#include <vtkPointData.h>
13#include <vtkPolyData.h>
14#include <vtkSmartPointer.h>
15
17
19 this->SetNumberOfInputPorts(0);
20 this->SetNumberOfOutputPorts(1);
21}
23
25 vtkInformation *ttkNotUsed(info)) {
26 return 0;
27}
28
29int ttkIcosphere::FillOutputPortInformation(int port, vtkInformation *info) {
30 if(port == 0)
31 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData");
32 else
33 return 0;
34 return 1;
35}
36
37int ttkIcosphere::RequestData(vtkInformation *ttkNotUsed(request),
38 vtkInformationVector **ttkNotUsed(inputVector),
39 vtkInformationVector *outputVector) {
40 // get parameter
41 size_t const nSpheres
42 = this->Centers ? this->Centers->GetNumberOfTuples() : 1;
43 size_t const nSubdivisions = this->GetNumberOfSubdivisions();
44 double const radius = this->GetRadius();
45
46 bool const useDoublePrecision
47 = this->Centers ? this->Centers->GetDataType() == VTK_DOUBLE : false;
48
49 // prepare the output buffers
50 size_t nVertices = 0;
51 size_t nTriangles = 0;
53 nVertices, nTriangles, nSubdivisions);
54
55 const size_t nTotalVertices = nSpheres * nVertices;
56 const size_t nTotalTriangles = nSpheres * nTriangles;
57
58 auto points = vtkSmartPointer<vtkPoints>::New();
59 points->SetDataType(useDoublePrecision ? VTK_DOUBLE : VTK_FLOAT);
60 points->SetNumberOfPoints(nTotalVertices);
61
63 if(this->ComputeNormals) {
64 if(useDoublePrecision)
66 else
68
69 normals->SetName("Normals");
70 normals->SetNumberOfComponents(3);
71 normals->SetNumberOfTuples(nTotalVertices);
72 }
73
75 offsets->SetNumberOfTuples(nTotalTriangles + 1);
76 auto offsetsData
77 = static_cast<vtkIdType *>(ttkUtils::GetVoidPointer(offsets));
78 for(size_t i = 0; i <= nTotalTriangles; i++)
79 offsetsData[i] = i * 3;
80
81 auto connectivity = vtkSmartPointer<vtkIdTypeArray>::New();
82 connectivity->SetNumberOfTuples(nTotalTriangles * 3);
83
84 int status = 0;
85 if(useDoublePrecision) {
86 using DT = double;
87 status = this->computeIcospheres<DT, vtkIdType>(
88 ttkUtils::GetPointer<DT>(points->GetData()),
89 ttkUtils::GetPointer<vtkIdType>(connectivity),
90
91 nSpheres, nSubdivisions, radius,
92 this->Centers ? ttkUtils::GetPointer<DT>(this->Centers) : this->Center,
93 this->ComputeNormals ? ttkUtils::GetPointer<DT>(normals) : nullptr);
94 } else {
95 using DT = float;
96 DT centerFloat[3]{
97 (DT)this->Center[0], (DT)this->Center[1], (DT)this->Center[2]};
98 status = this->computeIcospheres<DT, vtkIdType>(
99 ttkUtils::GetPointer<DT>(points->GetData()),
100 ttkUtils::GetPointer<vtkIdType>(connectivity),
101
102 nSpheres, nSubdivisions, radius,
103 this->Centers ? ttkUtils::GetPointer<DT>(this->Centers) : centerFloat,
104 this->ComputeNormals ? ttkUtils::GetPointer<DT>(normals) : nullptr);
105 }
106
107 if(!status)
108 return 0;
109
110 // finalize output
111 {
112 auto output = vtkPolyData::GetData(outputVector);
113 output->SetPoints(points);
114
116 cells->SetData(offsets, connectivity);
117 output->SetPolys(cells);
118
119 if(this->ComputeNormals)
120 output->GetPointData()->SetNormals(normals);
121 }
122
123 return 1;
124}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
~ttkIcosphere() override
int FillInputPortInformation(int port, vtkInformation *info) override
virtual double GetRadius()
virtual int GetNumberOfSubdivisions()
int FillOutputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int computeNumberOfVerticesAndTriangles(size_t &nVertices, size_t &nTriangles, const size_t nSubdivisions) const
Definition Icosphere.h:33
vtkStandardNewMacro(ttkIcosphere)