Flow123d  JS_before_hm-1003-g4e68d2c
output_vtk.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file output_vtk.cc
15  * @brief The functions for outputs to VTK files.
16  */
17 
18 #include "output_vtk.hh"
20 #include "element_data_cache.hh"
21 #include "output_mesh.hh"
22 #include "mesh/mesh.h"
23 
24 #include <limits.h>
25 #include "input/factory.hh"
27 #include "system/file_path.hh"
28 #include "tools/unit_si.hh"
29 #include "la/distribution.hh"
30 
31 #include "config.h"
32 
34 
35 
36 using namespace Input::Type;
37 
38 const Record & OutputVTK::get_input_type() {
39  return Record("vtk", "Parameters of vtk output format.")
40  // It is derived from abstract class
42  .declare_key("variant", OutputVTK::get_input_type_variant(), Default("\"ascii\""),
43  "Variant of output stream file format.")
44  // The parallel or serial variant
45  .declare_key("parallel", Bool(), Default("false"),
46  "Parallel or serial version of file format.")
47  .close();
48 }
49 
50 
52  return Selection("VTK variant (ascii or binary)")
54  "ASCII variant of VTK file format")
56  "Uncompressed appended binary XML VTK format without usage of base64 encoding of appended data.")
57 #ifdef FLOW123D_HAVE_ZLIB
59  "Appended binary XML VTK format without usage of base64 encoding of appended data. Compressed with ZLib.")
60 #endif // FLOW123D_HAVE_ZLIB
61  .close();
62 }
63 
64 
65 const int OutputVTK::registrar = Input::register_class< OutputVTK >("vtk") +
67 
68 
69 const std::vector<std::string> OutputVTK::formats = { "ascii", "appended", "appended" };
70 
71 
72 
74 {
75  this->enable_refinement_ = true;
76 }
77 
78 
79 
81 {
82  this->write_tail();
83 }
84 
85 
86 
87 void OutputVTK::init_from_input(const std::string &equation_name, const Input::Record &in_rec, std::string unit_str)
88 {
89  OutputTime::init_from_input(equation_name, in_rec, unit_str);
90 
91  auto format_rec = (Input::Record)(input_record_.val<Input::AbstractRecord>("format"));
92  variant_type_ = format_rec.val<VTKVariant>("variant");
93  this->parallel_ = format_rec.val<bool>("parallel");
94  this->fix_main_file_extension(".pvd");
95 
96  if(this->rank_ == 0) {
97  try {
98  this->_base_filename.open_stream( this->_base_file );
99  this->set_stream_precision(this->_base_file);
100  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, input_record_)
101 
102  LogOut() << "Writing flow output file: " << this->_base_filename << " ... ";
103  }
104 
105  this->make_subdirectory();
106  this->write_head();
107 
108 }
109 
110 string OutputVTK::form_vtu_filename_(string basename, int i_step, int rank) {
111  ostringstream ss;
112  if (this->parallel_) {
113  // parallel file
114  ss << basename << "/" << basename << "-"
115  << std::setw(6) << std::setfill('0') << i_step << "." << rank << ".vtu";
116  } else {
117  // serial file
118  ss << basename << "/" << basename << "-"
119  << std::setw(6) << std::setfill('0') << i_step << ".vtu";
120  }
121  return ss.str();
122 }
123 
124 string pvd_dataset_line(double step, int rank, string file) {
125  ostringstream ss;
126  ss
127  << "<DataSet timestep=\"" << step
128  << "\" group=\"\" part=\"" << rank
129  << "\" file=\"" << file
130  << "\"/>\n";
131  return ss.str();
132 }
133 
135 {
136  ASSERT_PTR(this->nodes_).error();
137 
138  /* Output of serial format is implemented only in the first process */
139  if ( (this->rank_ != 0) && (!parallel_) ) {
140  return 0;
141  }
142 
143  /**
144  * TODO:
145  * - common creation of the VTU filename
146  * - names of parallel file: <base name>_<frame>.<rank>.vtu
147  * - for serial case use rank=0
148  */
149 
150  /* Write DataSets to the PVD file only in the first process */
151  if (this->rank_ == 0) {
152  ASSERT(this->_base_file.is_open())(this->_base_filename).error();
153 
154  /* Set floating point precision to max */
155  //this->_base_file.precision(std::numeric_limits<double>::digits10);
156  //int current_step = this->get_parallel_current_step();
157 
158  /* Write dataset lines to the PVD file. */
159  double corrected_time = (isfinite(this->time)?this->time:0);
160  corrected_time /= UnitSI().s().convert_unit_from(this->unit_string_);
161  if (parallel_) {
162  for (int i_rank=0; i_rank<n_proc_; ++i_rank) {
163  string file = this->form_vtu_filename_(main_output_basename_, current_step, i_rank);
164  this->_base_file << pvd_dataset_line(corrected_time, i_rank, file);
165  }
166  } else {
167  string file = this->form_vtu_filename_(main_output_basename_, current_step, -1);
168  this->_base_file << pvd_dataset_line(corrected_time, 0, file);
169  }
170  }
171 
172  /* write VTU file */
173  {
174  /* Open VTU file */
175  std::string frame_file_name = this->form_vtu_filename_(main_output_basename_, current_step, this->rank_);
176  FilePath frame_file_path({main_output_dir_, frame_file_name}, FilePath::output_file);
177  try {
178  frame_file_path.open_stream(_data_file);
180  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, input_record_)
181 
182  LogOut() << __func__ << ": Writing output (frame: " << this->current_step
183  << ", rank: " << this->rank_
184  << ") file: " << frame_file_name << " ... ";
185 
186  this->write_vtk_vtu();
187 
188  /* Close stream for file of current frame */
189  _data_file.close();
190  //delete data_file;
191  //this->_data_file = NULL;
192 
193  LogOut() << "O.K.";
194 
195  }
196 
197  return 1;
198 }
199 
200 
201 
202 
204 {
205  ASSERT_EQ(this->_base_filename.extension(), ".pvd").error();
208 
209  if(this->rank_ == 0) {
210  vector<string> sub_path = { main_output_dir_, main_output_basename_, "__tmp__" };
211  FilePath fp(sub_path, FilePath::output_file);
212  fp.create_output_dir();
213  }
214 }
215 
216 
217 
218 
220 {
221  ofstream &file = this->_data_file;
222 
223  file << "<?xml version=\"1.0\"?>" << endl;
224  // TODO: test endianess of platform (this would be important, when raw
225  // data will be saved to the VTK file)
226  file << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\"";
227  if ( this->variant_type_ != VTKVariant::VARIANT_ASCII ) {
228  file << " header_type=\"UInt64\"";
229  }
230  if ( this->variant_type_ == VTKVariant::VARIANT_BINARY_ZLIB ) {
231  file << " compressor=\"vtkZLibDataCompressor\"";
232  }
233  file << ">" << endl;
234  file << "<UnstructuredGrid>" << endl;
235 }
236 
237 
238 
239 std::shared_ptr<ElementDataCache<unsigned int>> OutputVTK::fill_element_types_data()
240 {
241  auto &offsets = *( this->offsets_->get_component_data(0).get() );
242  unsigned int n_elements = offsets.size();
243 
244  auto types = std::make_shared<ElementDataCache<unsigned int>>("types", (unsigned int)ElementDataCacheBase::N_SCALAR, n_elements);
245  std::vector< unsigned int >& data = *( types->get_component_data(0).get() );
246  int n_nodes;
247 
248  n_nodes = offsets[0];
249  switch(n_nodes) {
250  case 2:
251  data[0] = (unsigned int)VTK_LINE;
252  break;
253  case 3:
254  data[0] = (unsigned int)VTK_TRIANGLE;
255  break;
256  case 4:
257  data[0] = (unsigned int)VTK_TETRA;
258  break;
259  }
260 
261  for(unsigned int i=1; i < n_elements; i++)
262  {
263  n_nodes = offsets[i]-offsets[i-1];
264  switch(n_nodes) {
265  case 2:
266  data[i] = (unsigned int)VTK_LINE;
267  break;
268  case 3:
269  data[i] = (unsigned int)VTK_TRIANGLE;
270  break;
271  case 4:
272  data[i] = (unsigned int)VTK_TETRA;
273  break;
274  }
275  }
276 
277  return types;
278 }
279 
280 
281 
283 {
284  // names of types in DataArray section
285  static const std::vector<std::string> types = {
286  "Int8", "UInt8", "Int16", "UInt16", "Int32", "UInt32", "Float32", "Float64" };
287 
288  ofstream &file = this->_data_file;
289 
290  file << "<DataArray type=\"" << types[output_data->vtk_type()] << "\" ";
291  // possibly write name
292  if( ! output_data->field_input_name().empty())
293  file << "Name=\"" << output_data->field_input_name() <<"\" ";
294  // write number of components
295  if (output_data->n_comp() > 1)
296  {
297  file
298  << "NumberOfComponents=\"" << output_data->n_comp() << "\" ";
299  }
300  file << "format=\"" << formats[this->variant_type_] << "\"";
301 
302  if ( this->variant_type_ == VTKVariant::VARIANT_ASCII ) {
303  // ascii output
304  file << ">" << endl;
305  //file << std::fixed << std::setprecision(10); // Set precision to max
306  output_data->print_ascii_all(file);
307  file << "\n</DataArray>" << endl;
308  } else {
309  // binary output is stored to appended_data_ stream
310  double range_min, range_max;
311  output_data->get_min_max_range(range_min, range_max);
312  file << " offset=\"" << appended_data_.tellp() << "\" ";
313  file << "RangeMin=\"" << range_min << "\" RangeMax=\"" << range_max << "\"/>" << endl;
314  if ( this->variant_type_ == VTKVariant::VARIANT_BINARY_UNCOMPRESSED ) {
315  output_data->print_binary_all( appended_data_ );
316  } else { // ZLib compression
317  stringstream uncompressed_data, compressed_data;
318  output_data->print_binary_all( uncompressed_data, false );
319  this->compress_data(uncompressed_data, compressed_data);
320  appended_data_ << compressed_data.str();
321  }
322  }
323 
324 }
325 
326 
327 void OutputVTK::compress_data(stringstream &uncompressed_stream, stringstream &compressed_stream) {
328  // size of block of compressed data.
329  static const size_t BUF_SIZE = 32 * 1024;
330 
331  string uncompressed_string = uncompressed_stream.str(); // full uncompressed string
332  zlib_ulong uncompressed_size = uncompressed_string.size(); // size of uncompressed string
333  stringstream compressed_data; // helper stream stores blocks of compress data
334 
335  zlib_ulong count_of_blocks = (uncompressed_size + BUF_SIZE - 1) / BUF_SIZE;
336  zlib_ulong last_block_size = (uncompressed_size % BUF_SIZE);
337  compressed_stream.write(reinterpret_cast<const char*>(&count_of_blocks), sizeof(unsigned long long int));
338  compressed_stream.write(reinterpret_cast<const char*>(&BUF_SIZE), sizeof(unsigned long long int));
339  compressed_stream.write(reinterpret_cast<const char*>(&last_block_size), sizeof(unsigned long long int));
340 
341  for (zlib_ulong i=0; i<count_of_blocks; ++i) {
342  // get block of data for compression
343  std::string data_block = uncompressed_string.substr(i*BUF_SIZE, BUF_SIZE);
344  zlib_ulong data_block_size = data_block.size();
345 
346  std::vector<uint8_t> buffer;
347  uint8_t temp_buffer[BUF_SIZE];
348 
349  // set zlib object
350  z_stream strm;
351  strm.zalloc = 0;
352  strm.zfree = 0;
353  strm.next_in = reinterpret_cast<uint8_t *>(&data_block[0]);
354  strm.avail_in = data_block_size;
355  strm.next_out = temp_buffer;
356  strm.avail_out = BUF_SIZE;
357 
358  // compression of data
359  deflateInit(&strm, Z_BEST_COMPRESSION);
360  while (strm.avail_in != 0) {
361  int res = deflate(&strm, Z_NO_FLUSH);
362  ASSERT_EQ(res, Z_OK).error();
363  if (strm.avail_out == 0) {
364  buffer.insert(buffer.end(), temp_buffer, temp_buffer + BUF_SIZE);
365  strm.next_out = temp_buffer;
366  strm.avail_out = BUF_SIZE;
367  }
368  }
369  int deflate_res = Z_OK;
370  while (deflate_res == Z_OK) {
371  if (strm.avail_out == 0) {
372  buffer.insert(buffer.end(), temp_buffer, temp_buffer + BUF_SIZE);
373  strm.next_out = temp_buffer;
374  strm.avail_out = BUF_SIZE;
375  }
376  deflate_res = deflate(&strm, Z_FINISH);
377  }
378  ASSERT_EQ(deflate_res, Z_STREAM_END).error();
379  buffer.insert(buffer.end(), temp_buffer, temp_buffer + BUF_SIZE - strm.avail_out);
380  deflateEnd(&strm);
381 
382  // store compress data and its size to streams
383  std::string str(buffer.begin(), buffer.end());
384  zlib_ulong compressed_data_size = str.size();
385  compressed_stream.write(reinterpret_cast<const char*>(&compressed_data_size), sizeof(unsigned long long int));
386  compressed_data << str;
387  }
388  // push compress data to returned stream
389  compressed_stream << compressed_data.str();
390 }
391 
392 
394 {
395  for(OutputDataPtr data : output_data_vec)
396  if( ! data->is_dummy())
397  write_vtk_data(data);
398 }
399 
400 
401 
402 
404  OutputDataFieldVec &output_data_vec)
405 {
406  if (output_data_vec.empty()) return;
407 
408  file << "Scalars=\"";
409  for(OutputDataPtr data : output_data_vec )
410  if (data->n_comp() == ElementDataCacheBase::N_SCALAR
411  && ! data->is_dummy())
412  file << data->field_input_name() << ",";
413  file << "\" ";
414 
415  file << "Vectors=\"";
416  for(OutputDataPtr data : output_data_vec )
417  if (data->n_comp() == ElementDataCacheBase::N_VECTOR
418  && ! data->is_dummy())
419  file << data->field_input_name() << ",";
420  file << "\" ";
421 
422  file << "Tensors=\"";
423  for(OutputDataPtr data : output_data_vec )
424  if (data->n_comp() == ElementDataCacheBase::N_TENSOR
425  && ! data->is_dummy())
426  file << data->field_input_name() << ",";
427  file << "\"";
428 }
429 
430 
432 {
433  ofstream &file = this->_data_file;
434 
435  // merge node and corner data
436  OutputDataFieldVec node_corner_data(output_data_vec_[NODE_DATA]);
437  node_corner_data.insert(node_corner_data.end(),
439 
440  if( ! node_corner_data.empty() ) {
441  /* Write <PointData begin */
442  file << "<PointData ";
443  write_vtk_data_names(file, node_corner_data);
444  file << ">" << endl;
445 
446  /* Write data on nodes */
447  this->write_vtk_field_data(output_data_vec_[NODE_DATA]);
448 
449  /* Write data in corners of elements */
451 
452  /* Write PointData end */
453  file << "</PointData>" << endl;
454  }
455 }
456 
457 
459 {
460  ofstream &file = this->_data_file;
461 
462  auto &data_map = this->output_data_vec_[ELEM_DATA];
463  if (data_map.empty()) return;
464 
465  /* Write CellData begin */
466  file << "<CellData ";
467  write_vtk_data_names(file, data_map);
468  file << ">" << endl;
469 
470  /* Write own data */
471  this->write_vtk_field_data(data_map);
472 
473  /* Write PointData end */
474  file << "</CellData>" << endl;
475 }
476 
477 
479 {
480  ofstream &file = this->_data_file;
481 
482  auto &data_map = this->output_data_vec_[NATIVE_DATA];
483  if (data_map.empty()) return;
484 
485  /* Write Flow123dData begin */
486  file << "<Flow123dData ";
487  write_vtk_data_names(file, data_map);
488  file << ">" << endl;
489 
490  /* Write own data */
491  for(OutputDataPtr output_data : data_map) {
492  file << "<DataArray type=\"Float64\" ";
493  file << "Name=\"" << output_data->field_input_name() <<"\" ";
494  file << "format=\"" << formats[this->variant_type_] << "\" ";
495  file << "dof_handler_hash=\"" << output_data->dof_handler_hash() << "\" ";
496  file << "n_dofs_per_element=\"" << output_data->n_comp() << "\"";
497 
498  if ( this->variant_type_ == VTKVariant::VARIANT_ASCII ) {
499  // ascii output
500  file << ">" << endl;
501  file << std::fixed << std::setprecision(10); // Set precision to max
502  output_data->print_ascii_all(file);
503  file << "\n</DataArray>" << endl;
504  } else {
505  // binary output is stored to appended_data_ stream
506  double range_min, range_max;
507  output_data->get_min_max_range(range_min, range_max);
508  file << " offset=\"" << appended_data_.tellp() << "\" ";
509  file << "RangeMin=\"" << range_min << "\" RangeMax=\"" << range_max << "\"/>" << endl;
510  if ( this->variant_type_ == VTKVariant::VARIANT_BINARY_UNCOMPRESSED ) {
511  output_data->print_binary_all( appended_data_ );
512  } else { // ZLib compression
513  stringstream uncompressed_data, compressed_data;
514  output_data->print_binary_all( uncompressed_data, false );
515  this->compress_data(uncompressed_data, compressed_data);
516  appended_data_ << compressed_data.str();
517  }
518  }
519  }
520 
521  /* Write Flow123dData end */
522  file << "</Flow123dData>" << endl;
523 }
524 
525 
527 {
528  ofstream &file = this->_data_file;
529 
530  file << "</UnstructuredGrid>" << endl;
531  if ( this->variant_type_ != VTKVariant::VARIANT_ASCII ) {
532  // appended data of binary compressed output
533  file << "<AppendedData encoding=\"raw\">" << endl;
534  // appended data starts with '_' character
535  file << "_" << appended_data_.str() << endl;
536  file << "</AppendedData>" << endl;
537  }
538  file << "</VTKFile>" << endl;
539 }
540 
541 
543 {
544  ofstream &file = this->_data_file;
545 
546  /* Write header */
547  this->write_vtk_vtu_head();
548 
549  /* Write Piece begin */
550  file << "<Piece NumberOfPoints=\"" << this->nodes_->n_values()
551  << "\" NumberOfCells=\"" << this->offsets_->n_values() <<"\">" << endl;
552 
553  /* Write VTK Geometry */
554  file << "<Points>" << endl;
555  write_vtk_data(this->nodes_);
556  file << "</Points>" << endl;
557 
558  /* Write VTK Topology */
559  file << "<Cells>" << endl;
561  write_vtk_data(this->offsets_);
562  auto types = fill_element_types_data();
563  write_vtk_data( types );
564  file << "</Cells>" << endl;
565 
566  /* Write VTK scalar and vector data on nodes to the file */
567  this->write_vtk_node_data();
568 
569  /* Write VTK data on elements */
570  this->write_vtk_element_data();
571 
572  /* Write own VTK native data (skipped by Paraview) */
573  this->write_vtk_native_data();
574 
575  /* Write Piece end */
576  file << "</Piece>" << endl;
577 
578  /* Write tail */
579  this->write_vtk_vtu_tail();
580 }
581 
582 
583 
585 {
586  /* Output to PVD file is implemented only in the first process */
587  if(this->rank_ != 0) {
588  return 0;
589  }
590 
591  LogOut() << __func__ << ": Writing output file (head) " << this->_base_filename << " ... ";
592 
593  this->_base_file << "<?xml version=\"1.0\"?>" << endl;
594  this->_base_file << "<VTKFile type=\"Collection\" version=\"0.1\" byte_order=\"LittleEndian\">" << endl;
595  this->_base_file << "<Collection>" << endl;
596 
597  LogOut() << "O.K.";
598 
599  return 1;
600 }
601 
602 
604 {
605  /* Output to PVD file is implemented only in the first process */
606  if(this->rank_ != 0) {
607  return 0;
608  }
609 
610  LogOut() << __func__ << ": Writing output file (tail) " << this->_base_filename << " ... ";
611 
612  this->_base_file << "</Collection>" << endl;
613  this->_base_file << "</VTKFile>" << endl;
614 
615  LogOut() << "O.K.";
616 
617  return 1;
618 }
619 
620 
621 
622 
623 
624 
Classes for auxiliary output mesh.
string stem() const
Definition: file_path.cc:193
double time
Definition: output_time.hh:275
Input::Record input_record_
Definition: output_time.hh:285
void fix_main_file_extension(std::string extension)
Definition: output_time.cc:150
uLong zlib_ulong
Definition: output_vtk.hh:47
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:602
void make_subdirectory()
Definition: output_vtk.cc:203
static const std::vector< std::string > formats
Formats of DataArray section.
Definition: output_vtk.hh:137
std::shared_ptr< ElementDataCache< unsigned int > > connectivity_
Vector maps the nodes to their coordinates in vector nodes_.
Definition: output_time.hh:325
static const Input::Type::Record & get_input_type()
The definition of input record for vtk file format.
Definition: output_vtk.cc:38
void create_output_dir()
Definition: file_path.cc:176
void write_vtk_field_data(OutputDataFieldVec &output_data_map)
Definition: output_vtk.cc:393
std::shared_ptr< ElementDataCache< unsigned int > > fill_element_types_data()
Fills the data cache with VTK element types indicators.
Definition: output_vtk.cc:239
void write_vtk_vtu(void)
This function write all scalar and vector data on nodes and elements to the VTK file (...
Definition: output_vtk.cc:542
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Class for declaration of the input of type Bool.
Definition: type_base.hh:452
std::shared_ptr< ElementDataCache< unsigned int > > offsets_
Vector of offsets of node indices of elements. Maps elements to their nodes in connectivity_.
Definition: output_time.hh:327
FilePath _base_filename
Definition: output_time.hh:295
void write_vtk_data(OutputDataPtr output_data)
Definition: output_vtk.cc:282
double convert_unit_from(std::string actual_unit) const
Convert and check user-defined unit.
Definition: unit_si.cc:217
Abstract linear system class.
Definition: balance.hh:40
static const int registrar
Registrar of class to factory.
Definition: output_vtk.hh:134
#define INPUT_CATCH(ExceptionType, AddressEITag, input_accessor)
Definition: accessors.hh:63
static Input::Type::Abstract & get_input_format_type()
The specification of output file format.
Definition: output_time.cc:64
std::shared_ptr< ElementDataCacheBase > OutputDataPtr
Definition: output_time.hh:122
std::shared_ptr< ElementDataCache< double > > nodes_
Vector of node coordinates. [spacedim x n_nodes].
Definition: output_time.hh:323
string main_output_basename_
Basename of main output file (without extension)
Definition: output_vtk.hh:229
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
int write_head(void)
This function writes header of VTK (.pvd) file format.
Definition: output_vtk.cc:584
string main_output_dir_
Main output file directory.
Definition: output_vtk.hh:232
void write_vtk_element_data(void)
Write data on elements to the VTK file (.vtu)
Definition: output_vtk.cc:458
#define LogOut()
Macro defining &#39;log&#39; record of log.
Definition: logger.hh:261
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
VTKVariant
The declaration enumeration used for variant of file VTK format.
Definition: output_vtk.hh:97
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
bool parallel_
Parallel or serial version of file format (parallel has effect only for VTK)
Definition: output_time.hh:317
int rank() const
Definition: output_time.hh:211
void open_stream(Stream &stream) const
Definition: file_path.cc:211
int write_data(void)
This function write data to VTK (.pvd) file format for curent time.
Definition: output_vtk.cc:134
ofstream _data_file
Definition: output_vtk.hh:216
int current_step
Definition: output_time.hh:270
ostringstream appended_data_
Definition: output_vtk.hh:221
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
const Ret val(const string &key) const
UnitSI & s(int exp=1)
Definition: unit_si.cc:76
void write_vtk_vtu_tail(void)
Write tail of VTK file (.vtu)
Definition: output_vtk.cc:526
This class is used for output data to VTK file format.
Definition: output_vtk.hh:43
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
string parent_path() const
Definition: file_path.cc:183
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:503
int write_tail(void)
This function writes tail of VTK (.pvd) file format.
Definition: output_vtk.cc:603
void set_stream_precision(std::ofstream &stream)
Definition: output_time.cc:99
string unit_string_
String representation of time unit.
Definition: output_time.hh:320
ofstream _base_file
Definition: output_time.hh:290
~OutputVTK()
The destructor of this class. It writes tail of the file too.
Definition: output_vtk.cc:80
string extension() const
Definition: file_path.cc:198
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:458
void write_vtk_native_data(void)
Write native data (part of our own data skipped by Paraview) to the VTK file (.vtu) ...
Definition: output_vtk.cc:478
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
virtual void init_from_input(const std::string &equation_name, const Input::Record &in_rec, std::string unit_str)
Constructor of OutputTime object. It opens base file for writing.
Definition: output_time.cc:83
Support classes for parallel programing.
OutputVTK()
The constructor of this class. The head of file is written, when constructor is called.
Definition: output_vtk.cc:73
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
const Selection & close() const
Close the Selection, no more values can be added.
string form_vtu_filename_(string basename, int i_step, int rank)
Definition: output_vtk.cc:110
VTKVariant variant_type_
Output format (ascii, binary or binary compressed)
Definition: output_vtk.hh:235
void write_vtk_node_data(void)
Write data on nodes to the VTK file (.vtu)
Definition: output_vtk.cc:431
string pvd_dataset_line(double step, int rank, string file)
Definition: output_vtk.cc:124
OutputDataFieldVec output_data_vec_[N_DISCRETE_SPACES]
Definition: output_time.hh:265
void compress_data(stringstream &uncompressed_stream, stringstream &compressed_stream)
Definition: output_vtk.cc:327
Record type proxy class.
Definition: type_record.hh:182
bool enable_refinement_
Auxiliary flag for refinement enabling, due to gmsh format.
Definition: output_time.hh:314
void write_vtk_vtu_head(void)
Write header of VTK file (.vtu)
Definition: output_vtk.cc:219
Class for representation SI units of Fields.
Definition: unit_si.hh:40
Template for classes storing finite set of named values.
void init_from_input(const std::string &equation_name, const Input::Record &in_rec, std::string unit_str) override
Override OutputTime::init_from_input.
Definition: output_vtk.cc:87
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:180
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual)
Definition: asserts.hh:328
static const Input::Type::Selection & get_input_type_variant()
The definition of input record for selection of variant of file format.
Definition: output_vtk.cc:51
void write_vtk_data_names(ofstream &file, OutputDataFieldVec &output_data_map)
Write names of data sets in output_data vector that have value type equal to type. Output is done into stream file.
Definition: output_vtk.cc:403