TTK
Loading...
Searching...
No Matches
ttkTopologicalCompressionWriter.cpp
Go to the documentation of this file.
1#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5#include <vtkDoubleArray.h>
6#include <vtkExecutive.h>
7#include <vtkFloatArray.h>
8#include <vtkIdTypeArray.h>
9#include <vtkImageData.h>
10#include <vtkInformation.h>
11#include <vtkIntArray.h>
12#include <vtkPointData.h>
13#include <vtkSignedCharArray.h>
14
16
18 SetNumberOfInputPorts(1);
19 this->setDebugMsgPrefix("TopologicalCompressionWriter");
20}
21
23 int port, vtkInformation *info) {
24 if(port == 0) {
25 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkImageData");
26 return 1;
27 }
28 return 0;
29}
30
32
33 this->printMsg("New writing task.");
34
35 if(ZFPOnly && ZFPTolerance < 0.0) {
36 this->printErr(
37 "Wrong ZFP absolute error tolerance for ZFP-only use, aborting.");
38 return 0;
39 }
40
41 vtkDataObject *input = GetInput();
42 vtkImageData *vti = vtkImageData::SafeDownCast(input);
43
44 auto inputScalarField = this->GetInputArrayToProcess(0, input);
45 if(inputScalarField == nullptr) {
46 return 0;
47 }
48
49 auto triangulation = ttkAlgorithm::GetTriangulation(vti);
50 if(triangulation == nullptr) {
51 return 0;
52 }
53 this->preconditionTriangulation(triangulation);
54
55 const auto inputOffsets
56 = this->GetOrderArray(vtkDataSet::SafeDownCast(input), 0, 1, false);
57
58 vtkSmartPointer<vtkDataArray> outputScalarField;
59
60 switch(inputScalarField->GetDataType()) {
61 case VTK_SIGNED_CHAR:
62 outputScalarField = vtkSmartPointer<vtkSignedCharArray>::New();
63 break;
64 case VTK_DOUBLE:
65 outputScalarField = vtkSmartPointer<vtkDoubleArray>::New();
66 break;
67 case VTK_FLOAT:
68 outputScalarField = vtkSmartPointer<vtkFloatArray>::New();
69 break;
70 case VTK_INT:
71 outputScalarField = vtkSmartPointer<vtkIntArray>::New();
72 break;
73 case VTK_ID_TYPE:
74 outputScalarField = vtkSmartPointer<vtkIdTypeArray>::New();
75 break;
76 default:
77 this->printErr("Unsupported data type :(");
78 // Do nothing.
79 return 0;
80 }
81
82 outputScalarField->SetNumberOfTuples(inputScalarField->GetNumberOfTuples());
83 outputScalarField->SetName(inputScalarField->GetName());
84 Modified();
85
86 // manage tolerance (relative % -> absolute)
87 std::array<double, 2> sfRange{};
88 inputScalarField->GetRange(sfRange.data());
89 this->relToAbsZFPTolerance(this->ZFPTolerance, sfRange);
90
92 inputScalarField->GetDataType(), triangulation->getType(),
93 this->execute(
94 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputScalarField)),
95 static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(inputOffsets)),
96 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(outputScalarField)),
97 *static_cast<TTK_TT *>(triangulation->getData())));
98
99 this->printMsg("Compression successful.");
100
101 // Open file.
102 FILE *fp;
103 if((fp = fopen(FileName, "wb")) == nullptr) {
104 this->printErr("System IO error while opening the file.");
105 return 0;
106 }
107
108 int dt = inputScalarField->GetDataType();
109 auto vp = static_cast<double *>(ttkUtils::GetVoidPointer(inputScalarField));
110
111 this->setFileName(FileName);
112 this->WriteToFile(fp, CompressionType, ZFPOnly, SQMethod.c_str(), dt,
113 vti->GetExtent(), vti->GetSpacing(), vti->GetOrigin(), vp,
114 Tolerance, ZFPTolerance, inputScalarField->GetName());
115
116 this->printMsg("Wrote to " + std::string{FileName} + ".");
117 return 1;
118}
119
121 // copied from ParaView's vtkWriter::GetInput()
122 if(this->GetNumberOfInputConnections(0) < 1) {
123 return nullptr;
124 }
125 return this->GetExecutive()->GetInputData(0, 0);
126}
127
129 // copied from ParaView's vtkWriter::SetInputData()
130 this->SetInputDataInternal(0, input);
131}
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
VTK-filter that wraps the topologicalCompressionWriter processing package.
int FillInputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition: ttkUtils.cpp:225
void setDebugMsgPrefix(const std::string &prefix)
Definition: Debug.h:364
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
Definition: Debug.h:118
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition: Debug.h:149
void relToAbsZFPTolerance(const double zfpRelTol, const std::array< double, 2 > &sfRange)
Switch from a relative (% of scalar field range) to an absolute ZFP tolerance.
int WriteToFile(FILE *fp, int compressionType, bool zfpOnly, const char *sqMethod, int dataType, int *dataExtent, double *dataSpacing, double *dataOrigin, double *data, double tolerance, double zfpTolerance, const std::string &dataArrayName)
void preconditionTriangulation(AbstractTriangulation *const triangulation)
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition: ttkMacros.h:69
vtkStandardNewMacro(ttkTopologicalCompressionWriter)