TTK
Loading...
Searching...
No Matches
ttkTriangulationReader.cpp
Go to the documentation of this file.
1#include <Triangulation.h>
3
4#include <vtkDataObject.h>
5#include <vtkInformation.h>
6#include <vtkObjectFactory.h>
7
8#include <ttkUtils.h>
9#include <vtkUnstructuredGrid.h>
10
12
14 this->setDebugMsgPrefix("TriangulationReader");
15
16 this->SetNumberOfInputPorts(1);
17 this->SetNumberOfOutputPorts(1);
18}
19
21 vtkInformation *info) {
22 if(port == 0) {
23 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkUnstructuredGrid");
24 return 1;
25 }
26 return 0;
27}
28
30 vtkInformation *info) {
31 if(port == 0) {
32 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
33 return 1;
34 }
35 return 0;
36}
37
39 if(this->TriangulationFilePath.length() < 8
40 || this->TriangulationFilePath
41 .substr(this->TriangulationFilePath.length() - 4, 4)
42 .compare(".tpt")
43 != 0) {
44 this->printErr(
45 "TTK Preconditioned Triangulation file has to end with '.tpt'.");
46 return 0;
47 }
48
49 return 1;
50}
51
53 vtkInformationVector **inputVector,
54 vtkInformationVector *outputVector) {
55 ttk::Timer timer;
56
57 auto input = vtkUnstructuredGrid::GetData(inputVector[0]);
58 auto output = vtkUnstructuredGrid::GetData(outputVector, 0);
59
60 if(input == nullptr || output == nullptr) {
61 this->printErr("Invalid input dataset");
62 return 0;
63 }
64
65 output->ShallowCopy(input);
66
67 auto triangulation = ttkAlgorithm::GetTriangulation(input);
68 if(triangulation == nullptr) {
69 this->printErr("Invalid input triangulation");
70 return 0;
71 }
72
73 auto explTri
74 = static_cast<ttk::ExplicitTriangulation *>(triangulation->getData());
75
76 if(!this->validateFilePath()) {
77 return 0;
78 }
79 std::ifstream in(this->TriangulationFilePath);
80 explTri->readFromFile(in);
81
82 this->printMsg("Restored triangulation from " + this->TriangulationFilePath,
83 1.0, timer.getElapsedTime(), 1);
84
85 return 1;
86}
#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)
TTK VTK-filter that reads a TTK Triangulation file.
int FillOutputPortInformation(int port, vtkInformation *info) override
int FillInputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
ExplicitTriangulation is a class that provides time efficient traversal methods on triangulations of ...
double getElapsedTime()
Definition Timer.h:15
vtkStandardNewMacro(ttkTriangulationReader)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)