TTK
Loading...
Searching...
No Matches
ttkCinemaWriter.cpp
Go to the documentation of this file.
1#include <ttkCinemaWriter.h>
2
3#include <vtkInformation.h>
4
5#include <vtkDataArray.h>
6#include <vtkDirectory.h>
7#include <vtkFieldData.h>
8#include <vtkImageData.h>
9#include <vtkObjectFactory.h>
10#include <vtkPointData.h>
11#include <vtkStdString.h>
12#include <vtkStringArray.h>
13#include <vtkTable.h>
14
15// writers common
16#include <vtkZLibDataCompressor.h>
17
18// CSV writers
19#include <vtkDelimitedTextReader.h>
20#include <vtkDelimitedTextWriter.h>
21#include <vtkMultiBlockDataSet.h>
22
23// product writers
24#include <vtkPNGWriter.h>
25#include <vtkXMLDataObjectWriter.h>
26#include <vtkXMLMultiBlockDataWriter.h>
27
28// file lock
29#include <boost/interprocess/sync/file_lock.hpp>
30#include <sys/stat.h>
31
33
35 this->setDebugMsgPrefix("CinemaWriter");
36
37 this->SetNumberOfInputPorts(1);
38 this->SetNumberOfOutputPorts(1);
39}
40
42
43int ttkCinemaWriter::FillInputPortInformation(int port, vtkInformation *info) {
44 if(port == 0) {
45 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataObject");
46 } else {
47 return 0;
48 }
49 return 1;
50}
51
52int ttkCinemaWriter::FillOutputPortInformation(int port, vtkInformation *info) {
53 if(port == 0) {
55 } else {
56 return 0;
57 }
58 return 1;
59}
60
61static int ensureFolder(const std::string &path) {
62 auto directory = vtkSmartPointer<vtkDirectory>::New();
63 if(directory->Open(path.data()) == 1
64 || vtkDirectory::MakeDirectory(path.data()) == 1)
65 return 1;
66 else
67 return 0;
68}
69
71 if(this->DatabasePath.length() < 4
72 || this->DatabasePath.substr(this->DatabasePath.length() - 4, 4)
73 .compare(".cdb")
74 != 0) {
75 this->printErr("Database path has to end with '.cdb'.");
76 return 0;
77 }
78
79 return 1;
80}
81
83 ttk::Timer t;
84 this->printMsg("Deleting CDB: " + this->DatabasePath, 0,
86
87 this->Modified();
88 if(this->ValidateDatabasePath() == 0)
89 return 0;
90 int const status = vtkDirectory::DeleteDirectory(this->DatabasePath.data());
91
92 this->printMsg("Deleting CDB: " + this->DatabasePath, 1, t.getElapsedTime());
93
94 return status;
95}
96
97int ttkCinemaWriter::GetLockFilePath(std::string &path) {
98 if(!this->ValidateDatabasePath())
99 return 0;
100
101 path = this->DatabasePath + ".lockfile";
102
103 return 1;
104}
105
107 std::string lockFilePath;
108 if(!this->GetLockFilePath(lockFilePath))
109 return 0;
110
111 std::ofstream output(lockFilePath);
112 output.close();
113
114 return 1;
115}
116
117// =============================================================================
118// Process Request
119// =============================================================================
120int ttkCinemaWriter::ProcessDataProduct(vtkDataObject *input) {
121
122 // ---------------------------------------------------------------------------
123 // Get Correct Data Product Extension
124 // ---------------------------------------------------------------------------
126 if(input->IsA("vtkDataSet")) {
128 vtkXMLDataObjectWriter::NewWriter(input->GetDataObjectType()));
129 } else if(input->IsA("vtkMultiBlockDataSet")) {
131 } else {
132 return 0;
133 }
134
135 xmlWriter->SetDataModeToAppended();
136 xmlWriter->SetCompressorTypeToZLib();
137 const auto compressor
138 = vtkZLibDataCompressor::SafeDownCast(xmlWriter->GetCompressor());
139 if(compressor != nullptr)
140 compressor->SetCompressionLevel(this->CompressionLevel);
141
142 std::string const productExtension = this->Format == FORMAT::VTK
143 ? xmlWriter->GetDefaultFileExtension()
144 : this->Format == FORMAT::PNG ? "png"
145 : "ttk";
146
147 // -------------------------------------------------------------------------
148 // Prepare Field Data
149 // -------------------------------------------------------------------------
150 auto inputFD = vtkSmartPointer<vtkFieldData>::New();
151 inputFD->ShallowCopy(input->GetFieldData());
152 size_t nFields = inputFD->GetNumberOfArrays();
153
154 // -------------------------------------------------------------------------
155 // Ignore Meta Fields
156 // -------------------------------------------------------------------------
157 {
158 std::vector<std::string> toIgnore;
159
160 // ignore file column
161 toIgnore.emplace_back("FILE");
162
163 // remove temporary columns
164 for(size_t i = 0; i < nFields; i++) {
165 std::string const name(inputFD->GetArrayName(i));
166 if(name.substr(0, 4).compare("_ttk") == 0)
167 toIgnore.emplace_back(name);
168 }
169
170 // delete columns from fd
171 for(const auto &name : toIgnore)
172 inputFD->RemoveArray(name.data());
173
174 nFields = inputFD->GetNumberOfArrays();
175 }
176
177 // ===========================================================================
178 // Determine ProductId and collect values
179 // ===========================================================================
180 std::string productId;
181 std::string rDataProductPath;
182 std::vector<std::string> fields;
183 std::vector<std::string> values;
184 {
185
186 if(nFields < 1) {
187 productId = "FILE";
188 } else {
189 for(size_t i = 0; i < nFields; i++) {
190 auto array = inputFD->GetAbstractArray(i);
191 auto n = array->GetNumberOfTuples();
192 auto m = array->GetNumberOfComponents();
193 std::string value;
194 if(n < 1) {
195 value = "null";
196 } else {
197 value = "";
198 for(int j = 0; j < n; j++) {
199 for(int k = 0; k < m; k++) {
200 value += array->GetVariantValue(j * m + k).ToString() + "|";
201 }
202 value = value.substr(0, value.size() - 1);
203 value += ";";
204 }
205 value = value.substr(0, value.size() - 1);
206 }
207
208 fields.emplace_back(array->GetName());
209 values.emplace_back(value);
210 }
211
212 productId = values[0];
213 for(size_t i = 1; i < nFields; i++)
214 productId += "_" + values[i];
215 }
216
217 rDataProductPath = "data/" + productId + "." + productExtension;
218 }
219
220 // print keys
221 {
222 std::vector<std::vector<std::string>> rows(nFields + 1);
223 for(size_t i = 0; i < nFields; i++)
224 rows[i] = {fields[i], values[i]};
225 rows[nFields] = {"Key", productId};
226
228 }
229
230 // ===========================================================================
231 // Update database
232 // ===========================================================================
233 {
234 // Initialize file lock for remaining operations
235 std::string lockFilePath;
236 if(!this->GetLockFilePath(lockFilePath))
237 return 0;
238
239 boost::interprocess::file_lock flock;
240 try {
241 flock = boost::interprocess::file_lock(lockFilePath.data());
242 flock.lock();
243 } catch(boost::interprocess::interprocess_exception &) {
244 }
245
246 std::string const csvPath = this->DatabasePath + "/data.csv";
247 struct stat info;
248
249 // -------------------------------------------------------------------------
250 // If data.csv file does not exist create it
251 // -------------------------------------------------------------------------
252 if(stat(csvPath.data(), &info) != 0) {
253 ttk::Timer t;
254 this->printMsg("Creating data.csv file", 0, ttk::debug::LineMode::REPLACE,
256
257 std::ofstream csvFile;
258 csvFile.open(csvPath.data());
259 if(!csvFile.is_open()) {
260 this->printErr("Unable to create 'data.csv' file.");
261 return 0;
262 }
263
264 std::string header;
265 std::string firstRow;
266 for(size_t i = 0; i < nFields; i++) {
267 header += fields[i] + ",";
268 firstRow += values[i] + ",";
269 }
270 header += "FILE";
271 firstRow += rDataProductPath;
272
273 csvFile << header << endl << firstRow << endl;
274
275 // Close file
276 csvFile.close();
277
278 this->printMsg("Creating data.csv file", 1, t.getElapsedTime(),
280 }
281
282 // -------------------------------------------------------------------------
283 // Update data.csv file
284 // -------------------------------------------------------------------------
285
286 // read data.csv file
287 auto csvTable = vtkSmartPointer<vtkTable>::New();
288 {
289 ttk::Timer t;
290 this->printMsg("Reading data.csv file", 0, ttk::debug::LineMode::REPLACE,
292
294 reader->SetFileName(csvPath.data());
295 reader->DetectNumericColumnsOff();
296 reader->SetHaveHeaders(true);
297 reader->SetFieldDelimiterCharacters(",");
298 reader->Update();
299 auto readerOutput = vtkTable::SafeDownCast(reader->GetOutput());
300 if(!readerOutput) {
301 this->printErr("Unable to read 'data.csv' file.");
302 return 0;
303 }
304
305 csvTable->ShallowCopy(readerOutput);
306
307 this->printMsg("Reading data.csv file", 1, t.getElapsedTime(),
309 }
310
311 // check CSV file integrity
312 std::vector<vtkStringArray *> fieldToCSVColumnMap(nFields);
313 size_t const nRows = csvTable->GetNumberOfRows();
314 size_t const nColumns = csvTable->GetNumberOfColumns();
315 {
316 // Check If CSV file is empty
317 if(nColumns == 0) {
318 this->printErr(
319 "Empty 'data.csv' file (vtkDelimitedTextReader limitation).");
320 return 0;
321 }
322
323 // Check If CSV file contains columns not present in fields
324 for(size_t i = 0; i < nColumns; i++) {
325 std::string const columnName = csvTable->GetColumnName(i);
326
327 // skip FILE column
328 if(columnName.compare("FILE") == 0)
329 continue;
330
331 bool exist = false;
332 for(size_t j = 0; j < nFields; j++) {
333 if(fields[j].compare(columnName) == 0)
334 exist = true;
335 }
336 if(!exist) {
337 this->printErr("'data.csv' file contains column '" + columnName
338 + "' not present in field data.");
339 this->printErr("Unable to insert data product into cinema database.");
340 return 0;
341 }
342 }
343
344 // check if data product has fields not present in the CSV file
345 for(size_t i = 0; i < nFields; i++) {
346 auto column = vtkStringArray::SafeDownCast(
347 csvTable->GetColumnByName(fields[i].data()));
348 if(!column) {
349 this->printErr("Data product has field data array '" + fields[i]
350 + "' no recorded in the data.csv file.");
351 return 0;
352 }
353 fieldToCSVColumnMap[i] = column;
354 }
355 }
356
357 // TODO: make dynamic
358 auto fileColumn
359 = vtkStringArray::SafeDownCast(csvTable->GetColumnByName("FILE"));
360 if(!fileColumn) {
361 this->printErr("'data.csv' file has no 'FILE' column");
362 return 0;
363 }
364
365 // -------------------------------------------------------------------------
366 // delete products of the database with same keys
367 // -------------------------------------------------------------------------
368 {
369 std::vector<int> rowsToDelete;
370 for(size_t i = 0; i < nRows; i++) {
371 auto equal = true;
372 for(size_t j = 0; j < nFields; j++)
373 if(values[j].compare(fieldToCSVColumnMap[j]->GetValue(i)) != 0)
374 equal = false;
375 if(equal)
376 rowsToDelete.emplace_back(i);
377 }
378
379 if(rowsToDelete.size() > 0) {
380 ttk::Timer t;
381 this->printMsg("Deleting products with same keys", 0,
384
385 for(int i = rowsToDelete.size() - 1; i >= 0; i--) {
386 auto path = fileColumn->GetValue(rowsToDelete[i]);
387
388 // Remove DataProduct
389 remove((this->DatabasePath + "/" + path).data());
390
391 // Remove Row from CSV
392 csvTable->RemoveRow(rowsToDelete[i]);
393 }
394
395 this->printMsg("Deleting products with same keys", 1,
398 }
399 }
400
401 // -----------------------------------------------------------------
402 // Update data.csv file
403 // -----------------------------------------------------------------
404 {
405 ttk::Timer t;
406 this->printMsg("Updating data.csv file", 0, ttk::debug::LineMode::REPLACE,
408
409 size_t const rowIndex = csvTable->GetNumberOfRows();
410 csvTable->InsertNextBlankRow();
411
412 for(size_t j = 0; j < nFields; j++)
413 fieldToCSVColumnMap[j]->SetValue(rowIndex, values[j]);
414
415 fileColumn->SetValue(rowIndex, rDataProductPath);
416
417 // Write data.csv file
419 csvWriter->SetUseStringDelimiter(false);
420 csvWriter->SetFileName((this->DatabasePath + "/data.csv").data());
421 csvWriter->SetInputData(csvTable);
422 csvWriter->Write();
423
424 this->printMsg("Updating data.csv file", 1, t.getElapsedTime(),
426 }
427 }
428
429 // =========================================================================
430 // Store Data products
431 // =========================================================================
432 {
433 // Write input to disk
434 ttk::Timer t;
435 this->printMsg("Writing data product to disk", 0,
437
438 switch(this->Format) {
439
440 case FORMAT::VTK: {
441 xmlWriter->SetFileName(
442 (this->DatabasePath + "/" + rDataProductPath).data());
443 xmlWriter->SetInputData(input);
444 xmlWriter->Write();
445 break;
446 }
447 case FORMAT::PNG: {
448 auto inputAsID = vtkImageData::SafeDownCast(input);
449 if(!inputAsID) {
450 this->printErr("PNG format requires input of type 'vtkImageData'.");
451 return 0;
452 }
453
454 // search color array
455 {
456 bool found = false;
457 auto inputPD = inputAsID->GetPointData();
458 for(int i = 0; i < inputPD->GetNumberOfArrays(); i++) {
459 auto array = inputPD->GetAbstractArray(i);
460 if(array->IsA("vtkUnsignedCharArray")) {
461 inputPD->SetActiveScalars(inputPD->GetArrayName(i));
462 found = true;
463 break;
464 }
465 }
466
467 if(!found) {
468 this->printErr("Input image does not have any color array.");
469 return 0;
470 }
471 }
472
473 auto imageWriter = vtkSmartPointer<vtkPNGWriter>::New();
474 imageWriter->SetCompressionLevel(this->CompressionLevel);
475 imageWriter->SetFileName(
476 (this->DatabasePath + "/" + rDataProductPath).data());
477 imageWriter->SetInputData(inputAsID);
478 imageWriter->Write();
479 break;
480 }
481 case FORMAT::TTK: {
482 // Topological Compression
483 if(!input->IsA("vtkImageData")) {
484 vtkErrorMacro(
485 "Cannot use Topological Compression without a vtkImageData");
486 return 0;
487 }
488
489 const auto inputData = vtkImageData::SafeDownCast(input);
490 const auto sf = this->GetInputArrayToProcess(0, inputData);
491
492 vtkNew<ttkTopologicalCompressionWriter> topologicalCompressionWriter{};
493 topologicalCompressionWriter->SetInputArrayToProcess(
494 0, 0, 0, 0, sf->GetName());
495
496 topologicalCompressionWriter->SetTolerance(this->Tolerance);
497 topologicalCompressionWriter->SetMaximumError(this->MaximumError);
498 topologicalCompressionWriter->SetZFPTolerance(this->ZFPTolerance);
499 topologicalCompressionWriter->SetCompressionType(this->CompressionType);
500 topologicalCompressionWriter->SetSQMethodPV(this->SQMethodPV);
501 topologicalCompressionWriter->SetZFPOnly(this->ZFPOnly);
502 topologicalCompressionWriter->SetSubdivide(this->Subdivide);
503 topologicalCompressionWriter->SetUseTopologicalSimplification(
504 this->UseTopologicalSimplification);
505
506 // Check that input scalar field is indeed scalar
507 if(sf->GetNumberOfComponents() != 1) {
508 vtkErrorMacro("Input scalar field should have only 1 component");
509 return 0;
510 }
511 topologicalCompressionWriter->SetDebugLevel(this->debugLevel_);
512 topologicalCompressionWriter->SetFileName(
513 (this->DatabasePath + "/" + rDataProductPath).data());
514 topologicalCompressionWriter->SetInputData(inputData);
515 topologicalCompressionWriter->Write();
516 break;
517 }
518 default:
519 this->printErr("Unsupported Format");
520 return 0;
521 }
522
523 this->printMsg("Writing data product to disk", 1, t.getElapsedTime(),
525 }
526 this->printMsg("Wrote " + productId + "." + productExtension);
528 return 1;
529}
530
531int ttkCinemaWriter::RequestData(vtkInformation *ttkNotUsed(request),
532 vtkInformationVector **inputVector,
533 vtkInformationVector *outputVector) {
534 ttk::Timer timer;
535
536 // Print Status
537 {
538 std::string format = this->Format == FORMAT::VTK ? "VTK"
539 : this->Format == FORMAT::PNG ? "PNG"
540 : "TTK";
541 this->printMsg({{"Database", this->DatabasePath},
542 {"C. Level", std::to_string(this->CompressionLevel)},
543 {"Format", format},
544 {"Iterate", this->IterateMultiBlock ? "Yes" : "No"}});
546 }
547
548 // -------------------------------------------------------------------------
549 // Copy Input to Output
550 // -------------------------------------------------------------------------
551 auto input = vtkDataObject::GetData(inputVector[0]);
552 auto output = vtkDataObject::GetData(outputVector);
553 if(this->ForwardInput)
554 output->ShallowCopy(input);
555
556 // -------------------------------------------------------------------------
557 // Prepare Database
558 // -------------------------------------------------------------------------
559 {
560 // Initialize file lock for remaining operations
561 std::string lockFilePath;
562 if(!this->GetLockFilePath(lockFilePath))
563 return 0;
564
565 boost::interprocess::file_lock flock;
566 try {
567 flock = boost::interprocess::file_lock(lockFilePath.data());
568 flock.lock();
569 } catch(boost::interprocess::interprocess_exception &) {
570 }
571
572 if(this->ValidateDatabasePath() == 0)
573 return 0;
574
575 if(ensureFolder(this->DatabasePath) == 0) {
576 this->printErr("Unable to open/create cinema database.");
577 return 0;
578 }
579
580 if(ensureFolder(this->DatabasePath + "/data") == 0) {
581 this->printErr("Unable to open/create cinema database.");
582 return 0;
583 }
584 }
585
586 auto inputAsMB = vtkMultiBlockDataSet::SafeDownCast(input);
587 if(this->IterateMultiBlock && inputAsMB) {
588 size_t const n = inputAsMB->GetNumberOfBlocks();
589 for(size_t i = 0; i < n; i++)
590 if(!this->ProcessDataProduct(inputAsMB->GetBlock(i)))
591 return 0;
592 } else if(!this->ProcessDataProduct(input))
593 return 0;
594
595 // Output Performance
596 {
597 std::string resultString = "Complete (#products: ";
598 resultString += !this->IterateMultiBlock || !inputAsMB
599 ? "1"
600 : std::to_string(inputAsMB->GetNumberOfBlocks());
601 resultString += ")";
602
603 // print stats
604 this->printMsg(resultString, 1, timer.getElapsedTime());
606 }
607
608 return 1;
609}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
TTK VTK-filter that writes input to disk.
int ProcessDataProduct(vtkDataObject *input)
~ttkCinemaWriter() override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int GetLockFilePath(std::string &path)
int debugLevel_
Definition Debug.h:379
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
double getElapsedTime()
Definition Timer.h:15
vtkStandardNewMacro(ttkCinemaWriter)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)