TTK
Loading...
Searching...
No Matches
ttkJacobiSet.cpp
Go to the documentation of this file.
1#include <vtkCellData.h>
2#include <vtkDataArray.h>
3#include <vtkDataSet.h>
4#include <vtkInformation.h>
5#include <vtkNew.h>
6#include <vtkObjectFactory.h>
7#include <vtkPointData.h>
8#include <vtkSignedCharArray.h>
9#include <vtkSmartPointer.h>
10#include <vtkUnstructuredGrid.h>
11
12#include <ttkJacobiSet.h>
13#include <ttkMacros.h>
14#include <ttkUtils.h>
15
16#include <array>
17
19
21 this->SetNumberOfInputPorts(1);
22 this->SetNumberOfOutputPorts(1);
23 ForceInputOffsetScalarField = false;
24}
25
26int ttkJacobiSet::FillInputPortInformation(int port, vtkInformation *info) {
27 if(port == 0) {
28 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
29 return 1;
30 }
31 return 0;
32}
33
34int ttkJacobiSet::FillOutputPortInformation(int port, vtkInformation *info) {
35 if(port == 0 || port == 1) {
36 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
37 return 1;
38 }
39 return 0;
40}
41
42template <class dataTypeU, class dataTypeV>
43int ttkJacobiSet::dispatch(const dataTypeU *const uField,
44 const dataTypeV *const vField,
45 ttk::Triangulation *const triangulation) {
47 triangulation->getType(),
48 this->execute(jacobiSet_, uField, vField,
49 *static_cast<TTK_TT *>(triangulation->getData()),
50 &isPareto_));
51 return 0;
52}
53
54int ttkJacobiSet::RequestData(vtkInformation *ttkNotUsed(request),
55 vtkInformationVector **inputVector,
56 vtkInformationVector *outputVector) {
57
58 const auto input = vtkDataSet::GetData(inputVector[0]);
59 auto output = vtkUnstructuredGrid::GetData(outputVector);
60
61 const auto uComponent = this->GetInputArrayToProcess(0, input);
62 const auto vComponent = this->GetInputArrayToProcess(1, input);
63
64 if(uComponent == nullptr || vComponent == nullptr)
65 return -1;
66
67 this->printMsg("U-component: `" + std::string{uComponent->GetName()} + "'");
68 this->printMsg("V-component: `" + std::string{vComponent->GetName()} + "'");
69
70 auto triangulation = ttkAlgorithm::GetTriangulation(input);
71 if(triangulation == nullptr)
72 return -1;
73 this->preconditionTriangulation(triangulation);
74
75 // point data
76 const auto offsetFieldU = this->GetOrderArray(
77 input, 0, triangulation, false, 2, ForceInputOffsetScalarField);
78 const auto offsetFieldV = this->GetOrderArray(
79 input, 1, triangulation, false, 3, ForceInputOffsetScalarField);
80
81 this->setSosOffsetsU(
82 static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(offsetFieldU)));
83 this->setSosOffsetsV(
84 static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(offsetFieldV)));
85
86#ifndef TTK_ENABLE_DOUBLE_TEMPLATING
87 if(uComponent->GetDataType() != vComponent->GetDataType()) {
88 this->printErr(
89 "Scalar fields should have same input type. Use TTKPointDataConverter or "
90 "TTKArrayEditor to convert array types.");
91 return 0;
92 }
93 switch(uComponent->GetDataType()) {
94 vtkTemplateMacro(
95 dispatch(static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(uComponent)),
96 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(vComponent)),
97 triangulation));
98 }
99#else
100 switch(vtkTemplate2PackMacro(
101 uComponent->GetDataType(), vComponent->GetDataType())) {
102 vtkTemplate2Macro(
103 dispatch(static_cast<VTK_T1 *>(ttkUtils::GetVoidPointer(uComponent)),
104 static_cast<VTK_T2 *>(ttkUtils::GetVoidPointer(vComponent)),
105 triangulation));
106 }
107#endif // TTK_ENABLE_DOUBLE_TEMPLATING
108
109 vtkNew<vtkSignedCharArray> edgeTypes{};
110
111 edgeTypes->SetNumberOfComponents(1);
112 edgeTypes->SetNumberOfTuples(2 * jacobiSet_.size());
113 edgeTypes->SetName("Critical Type");
114
115 vtkNew<vtkSignedCharArray> isPareto{};
116 isPareto->SetNumberOfComponents(1);
117 isPareto->SetNumberOfTuples(2 * jacobiSet_.size());
118 isPareto->SetName("IsPareto");
119
120 vtkNew<vtkPoints> pointSet{};
121 pointSet->SetNumberOfPoints(2 * jacobiSet_.size());
122
123 vtkNew<vtkCellArray> cellArray{};
124 vtkNew<vtkIdList> idList{};
125 idList->SetNumberOfIds(2);
126
127 size_t pointCount = 0;
128 std::array<double, 3> p{};
129 for(size_t i = 0; i < jacobiSet_.size(); i++) {
130
131 int const edgeId = jacobiSet_[i].first;
132 ttk::SimplexId vertexId0 = -1, vertexId1 = -1;
133 triangulation->getEdgeVertex(edgeId, 0, vertexId0);
134 triangulation->getEdgeVertex(edgeId, 1, vertexId1);
135
136 input->GetPoint(vertexId0, p.data());
137 pointSet->SetPoint(pointCount, p.data());
138 edgeTypes->SetTuple1(pointCount, (float)jacobiSet_[i].second);
139 isPareto->SetTuple1(pointCount, (float)isPareto_[i]);
140
141 idList->SetId(0, pointCount);
142 pointCount++;
143
144 input->GetPoint(vertexId1, p.data());
145 pointSet->SetPoint(pointCount, p.data());
146 edgeTypes->SetTuple1(pointCount, (float)jacobiSet_[i].second);
147 isPareto->SetTuple1(pointCount, (float)isPareto_[i]);
148 idList->SetId(1, pointCount);
149 pointCount++;
150
151 cellArray->InsertNextCell(idList);
152 }
153 output->SetPoints(pointSet);
154 output->SetCells(VTK_LINE, cellArray);
155 output->GetPointData()->AddArray(edgeTypes);
156 output->GetPointData()->AddArray(isPareto);
157
158 if(EdgeIds) {
159 vtkNew<ttkSimplexIdTypeArray> edgeIdArray{};
160 edgeIdArray->SetNumberOfComponents(1);
161 edgeIdArray->SetNumberOfTuples(jacobiSet_.size());
162 edgeIdArray->SetName("EdgeIds");
163
164 pointCount = 0;
165 for(size_t i = 0; i < jacobiSet_.size(); i++) {
166 edgeIdArray->SetTuple1(pointCount, (float)jacobiSet_[i].first);
167 pointCount++;
168 }
169
170 output->GetCellData()->AddArray(edgeIdArray);
171 } else {
172 output->GetCellData()->RemoveArray("EdgeIds");
173 }
174
175 if(VertexScalars) {
176
177 for(int i = 0; i < input->GetPointData()->GetNumberOfArrays(); i++) {
178
179 const auto scalarField = input->GetPointData()->GetArray(i);
180 vtkSmartPointer<vtkDataArray> const scalarArray{
181 scalarField->NewInstance()};
182
183 scalarArray->SetNumberOfComponents(scalarField->GetNumberOfComponents());
184 scalarArray->SetNumberOfTuples(2 * jacobiSet_.size());
185 scalarArray->SetName(scalarField->GetName());
186 std::vector<double> value(scalarField->GetNumberOfComponents());
187
188 for(size_t j = 0; j < jacobiSet_.size(); j++) {
189 int const edgeId = jacobiSet_[j].first;
190 ttk::SimplexId vertexId0 = -1, vertexId1 = -1;
191 triangulation->getEdgeVertex(edgeId, 0, vertexId0);
192 triangulation->getEdgeVertex(edgeId, 1, vertexId1);
193
194 scalarField->GetTuple(vertexId0, value.data());
195 scalarArray->SetTuple(2 * j, value.data());
196
197 scalarField->GetTuple(vertexId1, value.data());
198 scalarArray->SetTuple(2 * j + 1, value.data());
199 }
200 output->GetPointData()->AddArray(scalarArray);
201 }
202 } else {
203 for(int i = 0; i < input->GetPointData()->GetNumberOfArrays(); i++) {
204 output->GetPointData()->RemoveArray(
205 input->GetPointData()->GetArray(i)->GetName());
206 }
207 }
208
209 return 1;
210}
#define ttkTemplateMacro(triangulationType, call)
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, ttk::Triangulation *triangulation, const bool getGlobalOrder=false, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
TTK VTK-filter that computes the Jacobi set of a bivariate volumetric data-set.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
int dispatch(const dataTypeU *const uField, const dataTypeV *const vField, ttk::Triangulation *const triangulation)
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void setSosOffsetsV(const SimplexId *const sosOffsets)
Definition JacobiSet.h:103
void preconditionTriangulation(AbstractTriangulation *const triangulation)
Definition JacobiSet.h:118
void setSosOffsetsU(const SimplexId *const sosOffsets)
Definition JacobiSet.h:91
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
Triangulation::Type getType() const
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkJacobiSet)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)