Flow123d  JS_before_hm-1575-ga41e096
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  // Perform output of last time step
83  this->write_time_frame();
84 
85  this->write_tail();
86 }
87 
88 
89 
90 void OutputVTK::init_from_input(const std::string &equation_name, const Input::Record &in_rec, std::string unit_str)
91 {
92  OutputTime::init_from_input(equation_name, in_rec, unit_str);
93 
94  auto format_rec = (Input::Record)(input_record_.val<Input::AbstractRecord>("format"));
95  variant_type_ = format_rec.val<VTKVariant>("variant");
96  this->parallel_ = format_rec.val<bool>("parallel");
97  this->fix_main_file_extension(".pvd");
98 
99  if(this->rank_ == 0) {
100  try {
101  this->_base_filename.open_stream( this->_base_file );
102  this->set_stream_precision(this->_base_file);
103  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, input_record_)
104 
105  LogOut() << "Writing flow output file: " << this->_base_filename << " ... ";
106  }
107 
108  this->make_subdirectory();
109  this->write_head();
110 
111 }
112 
113 string OutputVTK::form_vtu_filename_(string basename, int i_step, int rank) {
114  ostringstream ss;
115  if (this->parallel_) {
116  // parallel file
117  ss << basename << "/" << basename << "-"
118  << std::setw(6) << std::setfill('0') << i_step << "." << rank << ".vtu";
119  } else {
120  // serial file
121  ss << basename << "/" << basename << "-"
122  << std::setw(6) << std::setfill('0') << i_step << ".vtu";
123  }
124  return ss.str();
125 }
126 
127 string pvd_dataset_line(double step, int rank, string file) {
128  ostringstream ss;
129  ss
130  << "<DataSet timestep=\"" << step
131  << "\" group=\"\" part=\"" << rank
132  << "\" file=\"" << file
133  << "\"/>\n";
134  return ss.str();
135 }
136 
138 {
139  ASSERT_PTR(this->nodes_).error();
140 
141  /* Output of serial format is implemented only in the first process */
142  if ( (this->rank_ != 0) && (!parallel_) ) {
143  return 0;
144  }
145 
146  /**
147  * TODO:
148  * - common creation of the VTU filename
149  * - names of parallel file: <base name>_<frame>.<rank>.vtu
150  * - for serial case use rank=0
151  */
152 
153  /* Write DataSets to the PVD file only in the first process */
154  if (this->rank_ == 0) {
155  ASSERT(this->_base_file.is_open())(this->_base_filename).error();
156 
157  /* Set floating point precision to max */
158  //this->_base_file.precision(std::numeric_limits<double>::digits10);
159  //int current_step = this->get_parallel_current_step();
160 
161  /* Write dataset lines to the PVD file. */
162  double corrected_time = (isfinite(this->registered_time_)?this->registered_time_:0);
163  corrected_time /= UnitSI().s().convert_unit_from(this->unit_string_);
164  if (parallel_) {
165  for (int i_rank=0; i_rank<n_proc_; ++i_rank) {
166  string file = this->form_vtu_filename_(main_output_basename_, current_step, i_rank);
167  this->_base_file << pvd_dataset_line(corrected_time, i_rank, file);
168  }
169  } else {
170  string file = this->form_vtu_filename_(main_output_basename_, current_step, -1);
171  this->_base_file << pvd_dataset_line(corrected_time, 0, file);
172  }
173  }
174 
175  /* write VTU file */
176  {
177  /* Open VTU file */
178  std::string frame_file_name = this->form_vtu_filename_(main_output_basename_, current_step, this->rank_);
179  FilePath frame_file_path({main_output_dir_, frame_file_name}, FilePath::output_file);
180  try {
181  frame_file_path.open_stream(_data_file);
183  } INPUT_CATCH(FilePath::ExcFileOpen, FilePath::EI_Address_String, input_record_)
184 
185  LogOut() << __func__ << ": Writing output (frame: " << this->current_step
186  << ", rank: " << this->rank_
187  << ") file: " << frame_file_name << " ... ";
188 
189  this->write_vtk_vtu();
190 
191  /* Close stream for file of current frame */
192  _data_file.close();
193  //delete data_file;
194  //this->_data_file = NULL;
195 
196  LogOut() << "O.K.";
197 
198  }
199 
200  return 1;
201 }
202 
203 
204 
205 
207 {
208  ASSERT_EQ(this->_base_filename.extension(), ".pvd").error();
211 
212  if(this->rank_ == 0) {
213  vector<string> sub_path = { main_output_dir_, main_output_basename_, "__tmp__" };
214  FilePath fp(sub_path, FilePath::output_file);
215  fp.create_output_dir();
216  }
217 }
218 
219 
220 
221 
223 {
224  ofstream &file = this->_data_file;
225 
226  file << "<?xml version=\"1.0\"?>" << endl;
227  // TODO: test endianess of platform (this would be important, when raw
228  // data will be saved to the VTK file)
229  file << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\"";
230  if ( this->variant_type_ != VTKVariant::VARIANT_ASCII ) {
231  file << " header_type=\"UInt64\"";
232  }
233  if ( this->variant_type_ == VTKVariant::VARIANT_BINARY_ZLIB ) {
234  file << " compressor=\"vtkZLibDataCompressor\"";
235  }
236  file << ">" << endl;
237  file << "<UnstructuredGrid>" << endl;
238 }
239 
240 
241 
242 std::shared_ptr<ElementDataCache<unsigned int>> OutputVTK::fill_element_types_data()
243 {
244  auto &offsets = *( this->offsets_->get_component_data(0).get() );
245  unsigned int n_elements = offsets.size()-1;
246 
247  auto types = std::make_shared<ElementDataCache<unsigned int>>("types", (unsigned int)ElementDataCacheBase::N_SCALAR, n_elements);
248  std::vector< unsigned int >& data = *( types->get_component_data(0).get() );
249  int n_nodes;
250 
251  for(unsigned int i=0; i < n_elements; i++)
252  {
253  n_nodes = offsets[i+1]-offsets[i];
254  switch(n_nodes) {
255  case 2:
256  data[i] = (unsigned int)VTK_LINE;
257  break;
258  case 3:
259  data[i] = (unsigned int)VTK_TRIANGLE;
260  break;
261  case 4:
262  data[i] = (unsigned int)VTK_TETRA;
263  break;
264  }
265  }
266 
267  return types;
268 }
269 
270 
271 
272 void OutputVTK::write_vtk_data(OutputTime::OutputDataPtr output_data, unsigned int start)
273 {
274  // names of types in DataArray section
275  static const std::vector<std::string> types = {
276  "Int8", "UInt8", "Int16", "UInt16", "Int32", "UInt32", "Float32", "Float64" };
277 
278  ofstream &file = this->_data_file;
279 
280  file << "<DataArray type=\"" << types[output_data->vtk_type()] << "\" ";
281  // possibly write name
282  if( ! output_data->field_input_name().empty())
283  file << "Name=\"" << output_data->field_input_name() <<"\" ";
284  // write number of components
285  if (output_data->n_comp() > 1)
286  {
287  file
288  << "NumberOfComponents=\"" << output_data->n_comp() << "\" ";
289  }
290  file << "format=\"" << formats[this->variant_type_] << "\"";
291 
292  if ( this->variant_type_ == VTKVariant::VARIANT_ASCII ) {
293  // ascii output
294  file << ">" << endl;
295  //file << std::fixed << std::setprecision(10); // Set precision to max
296  output_data->print_ascii_all(file, start);
297  file << "\n</DataArray>" << endl;
298  } else {
299  // binary output is stored to appended_data_ stream
300  double range_min, range_max;
301  output_data->get_min_max_range(range_min, range_max);
302  file << " offset=\"" << appended_data_.tellp() << "\" ";
303  file << "RangeMin=\"" << range_min << "\" RangeMax=\"" << range_max << "\"/>" << endl;
304  if ( this->variant_type_ == VTKVariant::VARIANT_BINARY_UNCOMPRESSED ) {
305  output_data->print_binary_all( appended_data_, true, start );
306  } else { // ZLib compression
307  stringstream uncompressed_data, compressed_data;
308  output_data->print_binary_all( uncompressed_data, false, start );
309  this->compress_data(uncompressed_data, compressed_data);
310  appended_data_ << compressed_data.str();
311  }
312  }
313 
314 }
315 
316 
317 void OutputVTK::compress_data(stringstream &uncompressed_stream, stringstream &compressed_stream) {
318  // size of block of compressed data.
319  static const size_t BUF_SIZE = 32 * 1024;
320 
321  string uncompressed_string = uncompressed_stream.str(); // full uncompressed string
322  zlib_ulong uncompressed_size = uncompressed_string.size(); // size of uncompressed string
323  stringstream compressed_data; // helper stream stores blocks of compress data
324 
325  zlib_ulong count_of_blocks = (uncompressed_size + BUF_SIZE - 1) / BUF_SIZE;
326  zlib_ulong last_block_size = (uncompressed_size % BUF_SIZE);
327  compressed_stream.write(reinterpret_cast<const char*>(&count_of_blocks), sizeof(unsigned long long int));
328  compressed_stream.write(reinterpret_cast<const char*>(&BUF_SIZE), sizeof(unsigned long long int));
329  compressed_stream.write(reinterpret_cast<const char*>(&last_block_size), sizeof(unsigned long long int));
330 
331  for (zlib_ulong i=0; i<count_of_blocks; ++i) {
332  // get block of data for compression
333  std::string data_block = uncompressed_string.substr(i*BUF_SIZE, BUF_SIZE);
334  zlib_ulong data_block_size = data_block.size();
335 
336  std::vector<uint8_t> buffer;
337  uint8_t temp_buffer[BUF_SIZE];
338 
339  // set zlib object
340  z_stream strm;
341  strm.zalloc = 0;
342  strm.zfree = 0;
343  strm.next_in = reinterpret_cast<uint8_t *>(&data_block[0]);
344  strm.avail_in = data_block_size;
345  strm.next_out = temp_buffer;
346  strm.avail_out = BUF_SIZE;
347 
348  // compression of data
349  deflateInit(&strm, Z_BEST_COMPRESSION);
350  while (strm.avail_in != 0) {
351  int res = deflate(&strm, Z_NO_FLUSH);
352  ASSERT_EQ(res, Z_OK).error();
353  if (strm.avail_out == 0) {
354  buffer.insert(buffer.end(), temp_buffer, temp_buffer + BUF_SIZE);
355  strm.next_out = temp_buffer;
356  strm.avail_out = BUF_SIZE;
357  }
358  }
359  int deflate_res = Z_OK;
360  while (deflate_res == Z_OK) {
361  if (strm.avail_out == 0) {
362  buffer.insert(buffer.end(), temp_buffer, temp_buffer + BUF_SIZE);
363  strm.next_out = temp_buffer;
364  strm.avail_out = BUF_SIZE;
365  }
366  deflate_res = deflate(&strm, Z_FINISH);
367  }
368  ASSERT_EQ(deflate_res, Z_STREAM_END).error();
369  buffer.insert(buffer.end(), temp_buffer, temp_buffer + BUF_SIZE - strm.avail_out);
370  deflateEnd(&strm);
371 
372  // store compress data and its size to streams
373  std::string str(buffer.begin(), buffer.end());
374  zlib_ulong compressed_data_size = str.size();
375  compressed_stream.write(reinterpret_cast<const char*>(&compressed_data_size), sizeof(unsigned long long int));
376  compressed_data << str;
377  }
378  // push compress data to returned stream
379  compressed_stream << compressed_data.str();
380 }
381 
382 
384 {
385  for(OutputDataPtr data : output_data_vec)
386  if( ! data->is_dummy())
387  write_vtk_data(data);
388 }
389 
390 
391 
392 
394  OutputDataFieldVec &output_data_vec)
395 {
396  if (output_data_vec.empty()) return;
397 
398  file << "Scalars=\"";
399  for(OutputDataPtr data : output_data_vec )
400  if (data->n_comp() == ElementDataCacheBase::N_SCALAR
401  && ! data->is_dummy())
402  file << data->field_input_name() << ",";
403  file << "\" ";
404 
405  file << "Vectors=\"";
406  for(OutputDataPtr data : output_data_vec )
407  if (data->n_comp() == ElementDataCacheBase::N_VECTOR
408  && ! data->is_dummy())
409  file << data->field_input_name() << ",";
410  file << "\" ";
411 
412  file << "Tensors=\"";
413  for(OutputDataPtr data : output_data_vec )
414  if (data->n_comp() == ElementDataCacheBase::N_TENSOR
415  && ! data->is_dummy())
416  file << data->field_input_name() << ",";
417  file << "\"";
418 }
419 
420 
422 {
423  ofstream &file = this->_data_file;
424 
425  // merge node and corner data
426  OutputDataFieldVec node_corner_data(output_data_vec_[NODE_DATA]);
427  node_corner_data.insert(node_corner_data.end(),
429 
430  if( ! node_corner_data.empty() ) {
431  /* Write <PointData begin */
432  file << "<PointData ";
433  write_vtk_data_names(file, node_corner_data);
434  file << ">" << endl;
435 
436  /* Write data on nodes */
437  this->write_vtk_field_data(output_data_vec_[NODE_DATA]);
438 
439  /* Write data in corners of elements */
441 
442  /* Write PointData end */
443  file << "</PointData>" << endl;
444  }
445 }
446 
447 
449 {
450  ofstream &file = this->_data_file;
451 
452  auto &data_map = this->output_data_vec_[ELEM_DATA];
453  if (data_map.empty()) return;
454 
455  /* Write CellData begin */
456  file << "<CellData ";
457  write_vtk_data_names(file, data_map);
458  file << ">" << endl;
459 
460  /* Write own data */
461  this->write_vtk_field_data(data_map);
462 
463  /* Write PointData end */
464  file << "</CellData>" << endl;
465 }
466 
467 
469 {
470  ofstream &file = this->_data_file;
471 
472  auto &data_map = this->output_data_vec_[NATIVE_DATA];
473  if (data_map.empty()) return;
474 
475  /* Write Flow123dData begin */
476  file << "<Flow123dData ";
477  write_vtk_data_names(file, data_map);
478  file << ">" << endl;
479 
480  /* Write own data */
481  for(OutputDataPtr output_data : data_map) {
482  file << "<DataArray type=\"Float64\" ";
483  file << "Name=\"" << output_data->field_input_name() <<"\" ";
484  file << "format=\"" << formats[this->variant_type_] << "\" ";
485  file << "dof_handler_hash=\"" << output_data->dof_handler_hash() << "\" ";
486  file << "n_dofs_per_element=\"" << output_data->n_comp() << "\"";
487 
488  if ( this->variant_type_ == VTKVariant::VARIANT_ASCII ) {
489  // ascii output
490  file << ">" << endl;
491  file << std::fixed << std::setprecision(10); // Set precision to max
492  output_data->print_ascii_all(file);
493  file << "\n</DataArray>" << endl;
494  } else {
495  // binary output is stored to appended_data_ stream
496  double range_min, range_max;
497  output_data->get_min_max_range(range_min, range_max);
498  file << " offset=\"" << appended_data_.tellp() << "\" ";
499  file << "RangeMin=\"" << range_min << "\" RangeMax=\"" << range_max << "\"/>" << endl;
500  if ( this->variant_type_ == VTKVariant::VARIANT_BINARY_UNCOMPRESSED ) {
501  output_data->print_binary_all( appended_data_ );
502  } else { // ZLib compression
503  stringstream uncompressed_data, compressed_data;
504  output_data->print_binary_all( uncompressed_data, false );
505  this->compress_data(uncompressed_data, compressed_data);
506  appended_data_ << compressed_data.str();
507  }
508  }
509  }
510 
511  /* Write Flow123dData end */
512  file << "</Flow123dData>" << endl;
513 }
514 
515 
517 {
518  ofstream &file = this->_data_file;
519 
520  file << "</UnstructuredGrid>" << endl;
521  if ( this->variant_type_ != VTKVariant::VARIANT_ASCII ) {
522  // appended data of binary compressed output
523  file << "<AppendedData encoding=\"raw\">" << endl;
524  // appended data starts with '_' character
525  file << "_" << appended_data_.str() << endl;
526  file << "</AppendedData>" << endl;
527  }
528  file << "</VTKFile>" << endl;
529 }
530 
531 
533 {
534  ofstream &file = this->_data_file;
535 
536  /* Write header */
537  this->write_vtk_vtu_head();
538 
539  /* Write Piece begin */
540  file << "<Piece NumberOfPoints=\"" << this->nodes_->n_values()
541  << "\" NumberOfCells=\"" << this->offsets_->n_values()-1 <<"\">" << endl;
542 
543  /* Write VTK Geometry */
544  file << "<Points>" << endl;
545  write_vtk_data(this->nodes_);
546  file << "</Points>" << endl;
547 
548  /* Write VTK Topology */
549  file << "<Cells>" << endl;
551  write_vtk_data(this->offsets_, 1);
552  auto types = fill_element_types_data();
553  write_vtk_data( types );
554  file << "</Cells>" << endl;
555 
556  /* Write VTK scalar and vector data on nodes to the file */
557  this->write_vtk_node_data();
558 
559  /* Write VTK data on elements */
560  this->write_vtk_element_data();
561 
562  /* Write own VTK native data (skipped by Paraview) */
563  this->write_vtk_native_data();
564 
565  /* Write Piece end */
566  file << "</Piece>" << endl;
567 
568  /* Write tail */
569  this->write_vtk_vtu_tail();
570 }
571 
572 
573 
575 {
576  /* Output to PVD file is implemented only in the first process */
577  if(this->rank_ != 0) {
578  return 0;
579  }
580 
581  LogOut() << __func__ << ": Writing output file (head) " << this->_base_filename << " ... ";
582 
583  this->_base_file << "<?xml version=\"1.0\"?>" << endl;
584  this->_base_file << "<VTKFile type=\"Collection\" version=\"0.1\" byte_order=\"LittleEndian\">" << endl;
585  this->_base_file << "<Collection>" << endl;
586 
587  LogOut() << "O.K.";
588 
589  return 1;
590 }
591 
592 
594 {
595  /* Output to PVD file is implemented only in the first process */
596  if(this->rank_ != 0) {
597  return 0;
598  }
599 
600  LogOut() << __func__ << ": Writing output file (tail) " << this->_base_filename << " ... ";
601 
602  this->_base_file << "</Collection>" << endl;
603  this->_base_file << "</VTKFile>" << endl;
604 
605  LogOut() << "O.K.";
606 
607  return 1;
608 }
609 
610 
611 
612 
613 
614 
Classes for auxiliary output mesh.
string stem() const
Definition: file_path.cc:193
Input::Record input_record_
Definition: output_time.hh:290
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:206
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:330
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:383
std::shared_ptr< ElementDataCache< unsigned int > > fill_element_types_data()
Fills the data cache with VTK element types indicators.
Definition: output_vtk.cc:242
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:532
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:332
FilePath _base_filename
Definition: output_time.hh:300
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
void write_vtk_data(OutputDataPtr output_data, unsigned int start=0)
Definition: output_vtk.cc:272
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:328
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:574
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:448
#define LogOut()
Macro defining &#39;log&#39; record of log.
Definition: logger.hh:273
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:322
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:137
ofstream _data_file
Definition: output_vtk.hh:216
int current_step
Definition: output_time.hh:275
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:516
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:593
void set_stream_precision(std::ofstream &stream)
Definition: output_time.cc:99
string unit_string_
String representation of time unit.
Definition: output_time.hh:325
ofstream _base_file
Definition: output_time.hh:295
~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:468
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
void write_time_frame()
Definition: output_time.cc:198
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:113
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:421
string pvd_dataset_line(double step, int rank, string file)
Definition: output_vtk.cc:127
OutputDataFieldVec output_data_vec_[N_DISCRETE_SPACES]
Definition: output_time.hh:270
void compress_data(stringstream &uncompressed_stream, stringstream &compressed_stream)
Definition: output_vtk.cc:317
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:319
void write_vtk_vtu_head(void)
Write header of VTK file (.vtu)
Definition: output_vtk.cc:222
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:90
#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
double registered_time_
Definition: output_time.hh:280
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:393