TTK
Loading...
Searching...
No Matches
ttkTopologicalCompressionReader.cpp
Go to the documentation of this file.
1#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5#include <vtkDataObject.h>
6#include <vtkDataSet.h>
7#include <vtkDoubleArray.h>
8#include <vtkIdTypeArray.h>
9#include <vtkImageData.h>
10#include <vtkInformation.h>
11#include <vtkInformationVector.h>
12#include <vtkIntArray.h>
13#include <vtkNew.h>
14#include <vtkPointData.h>
15#include <vtkStreamingDemandDrivenPipeline.h>
16
18
20 SetNumberOfInputPorts(0);
21 SetNumberOfOutputPorts(1);
22 this->setDebugMsgPrefix("TopologicalCompressionReader");
23}
24
26 int port, vtkInformation *info) {
27 if(port == 0) {
28 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkImageData");
29 return 1;
30 }
31 return 0;
32}
33
35 vtkInformation *ttkNotUsed(request),
36 vtkInformationVector **ttkNotUsed(inputVector),
37 vtkInformationVector *outputVector) {
38
39 if(FileName == nullptr) {
40 return 1;
41 }
42
43 FILE *fp = fopen(FileName, "rb"); // binary mode
44 if(fp == nullptr) {
45 return 1;
46 }
47
48 // Fill spacing, origin, extent, scalar type
49 // L8 tolerance, ZFP factor
50 const auto res = this->ReadMetaData(fp);
51 if(res != 0) {
52 return 1;
53 }
54 DataScalarType = this->getDataScalarType();
55 for(int i = 0; i < 3; ++i) {
56 DataSpacing[i] = this->getDataSpacing()[i];
57 DataOrigin[i] = this->getDataOrigin()[i];
58 DataExtent[i] = this->getDataExtent()[i];
59 DataExtent[3 + i] = this->getDataExtent()[3 + i];
60 }
61
62 // ReadMetaData(fp);
63
64 vtkInformation *outInfo = outputVector->GetInformationObject(0);
65 outInfo->Set(vtkDataObject::SPACING(), DataSpacing.data(), 3);
66 outInfo->Set(vtkDataObject::ORIGIN(), DataOrigin.data(), 3);
67 outInfo->Set(
68 vtkStreamingDemandDrivenPipeline::WHOLE_EXTENT(), DataExtent.data(), 6);
69
70 int numberOfVertices = 1;
71 for(int i = 0; i < 3; ++i)
72 numberOfVertices *= (1 + DataExtent[2 * i + 1] - DataExtent[2 * i]);
73 outInfo->Set(vtkDataObject::FIELD_NUMBER_OF_TUPLES(), numberOfVertices);
74
75 vtkDataObject::SetPointDataActiveScalarInfo(outInfo, DataScalarType, 1);
76
77 rewind(fp);
78 fclose(fp);
79 fp = nullptr;
80
81 return 1;
82}
83
85 vtkInformation *ttkNotUsed(request),
86 vtkInformationVector **ttkNotUsed(inputVector),
87 vtkInformationVector *outputVector) {
88
89 // Initialize
90 if(FileName == nullptr) {
91 return 1;
92 }
93 FILE *fp = fopen(FileName, "rb"); // binary mode
94
95 if(fp == nullptr) {
96 return 1;
97 }
98
99 this->setFileName(FileName);
100 const auto res = this->ReadMetaData(fp);
101 if(res != 0) {
102 return 1;
103 }
104 DataScalarType = this->getDataScalarType();
105 for(int i = 0; i < 3; ++i) {
106 DataSpacing[i] = this->getDataSpacing()[i];
107 DataOrigin[i] = this->getDataOrigin()[i];
108 DataExtent[i] = this->getDataExtent()[i];
109 DataExtent[3 + i] = this->getDataExtent()[3 + i];
110 }
111 int const nx = 1 + DataExtent[1] - DataExtent[0];
112 int const ny = 1 + DataExtent[3] - DataExtent[2];
113 int const nz = 1 + DataExtent[5] - DataExtent[4];
114 int const vertexNumber = nx * ny * nz;
115 ZFPOnly = this->getZFPOnly();
116
117 vtkNew<vtkImageData> mesh{};
118 BuildMesh(mesh);
119
120 auto triangulation = ttkAlgorithm::GetTriangulation(mesh);
121 this->preconditionTriangulation(triangulation);
122
123 int status{0};
124 ttkTemplateMacro(triangulation->getType(),
125 status = this->ReadFromFile(
126 fp, *static_cast<TTK_TT *>(triangulation->getData())));
127 if(status != 0) {
128 vtkWarningMacro("Failure when reading compressed TTK file");
129 }
130
131 mesh->GetPointData()->RemoveArray(0);
132 mesh->GetPointData()->SetNumberOfTuples(vertexNumber);
133
134 vtkNew<vtkDoubleArray> decompressed{};
135 decompressed->SetNumberOfTuples(vertexNumber);
136 const auto &name = this->getDataArrayName();
137 if(!name.empty()) {
138 decompressed->SetName(name.data());
139 } else {
140 decompressed->SetName("Decompressed");
141 }
142 for(int i = 0; i < vertexNumber; ++i) {
143 decompressed->SetTuple1(i, decompressedData_[i]);
144 }
145
146 // decompressed->SetVoidArray(, vertexNumber, 0);
147 mesh->GetPointData()->AddArray(decompressed);
148
149 if(SQMethodInt != 1 && SQMethodInt != 2 && !ZFPOnly) {
150 vtkNew<vtkIntArray> vertexOffset{};
151 vertexOffset->SetNumberOfTuples(vertexNumber);
152 vertexOffset->SetName(this->GetOrderArrayName(decompressed).data());
153 const auto &voidOffsets = this->getDecompressedOffsets();
154 for(size_t i = 0; i < voidOffsets.size(); ++i)
155 vertexOffset->SetTuple1(i, voidOffsets[i]);
156 // vertexOffset->SetVoidArray(
157 // topologicalCompression.getDecompressedOffsets(), vertexNumber, 0);
158 mesh->GetPointData()->AddArray(vertexOffset);
159 }
160
161 this->printMsg("Read " + std::to_string(mesh->GetNumberOfPoints())
162 + " vertice(s), " + std::to_string(mesh->GetNumberOfCells())
163 + " cell(s).");
164
165 // get the info object
166 vtkInformation *outInfo = outputVector->GetInformationObject(0);
167 outInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_NUMBER_OF_PIECES(), 1);
168
169 // Set the output
170 auto output = vtkImageData::GetData(outputVector);
171 output->ShallowCopy(mesh);
172
173 return 1;
174}
175
177 // copied from ParaView's vtkImageAlgorithm::GetOutput(int port)
178 return vtkImageData::SafeDownCast(this->GetOutputDataObject(0));
179}
180
181void ttkTopologicalCompressionReader::BuildMesh(vtkImageData *mesh) const {
182 int const nx = 1 + DataExtent[1] - DataExtent[0];
183 int const ny = 1 + DataExtent[3] - DataExtent[2];
184 int const nz = 1 + DataExtent[5] - DataExtent[4];
185 mesh->SetDimensions(nx, ny, nz);
186 mesh->SetSpacing(DataSpacing[0], DataSpacing[1], DataSpacing[2]);
187 mesh->SetOrigin(DataOrigin[0], DataOrigin[1], DataOrigin[2]);
188 mesh->AllocateScalars(DataScalarType, 2);
189 mesh->GetPointData()->SetNumberOfTuples(nx * ny * nz);
190}
#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)
static std::string GetOrderArrayName(vtkDataArray *const array)
VTK-filter that wraps the topologicalCompressionWriter processing package.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int, vtkInformation *) override
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
const std::vector< char > & getDataArrayName() const
std::vector< double > decompressedData_
std::vector< SimplexId > & getDecompressedOffsets()
void preconditionTriangulation(AbstractTriangulation *const triangulation)
std::string to_string(__int128)
Definition ripserpy.cpp:99
vtkStandardNewMacro(ttkTopologicalCompressionReader)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)