Flow123d  last_with_con_2.0.0-4-g42e6930
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"
19 #include "output_data_base.hh"
20 #include "output_mesh_data.hh"
21 #include "output_mesh.hh"
22 
23 #include <limits.h>
24 #include "input/factory.hh"
26 #include "system/file_path.hh"
27 
29 
30 
31 using namespace Input::Type;
32 
33 const Record & OutputVTK::get_input_type() {
34  return Record("vtk", "Parameters of vtk output format.")
35  // It is derived from abstract class
37  .declare_key("variant", OutputVTK::get_input_type_variant(), Default("\"ascii\""),
38  "Variant of output stream file format.")
39  // The parallel or serial variant
40  .declare_key("parallel", Bool(), Default("false"),
41  "Parallel or serial version of file format.")
42  // Type of compression
43  .declare_key("compression", OutputVTK::get_input_type_compression(), Default("\"none\""),
44  "Compression used in output stream file format.")
45  .close();
46 }
47 
48 
50  return Selection("VTK variant (ascii or binary)")
52  "ASCII variant of VTK file format")
54  "Binary variant of VTK file format (not supported yet)")
55  .close();
56 }
57 
58 
60  return Selection("Type of compression of VTK file format")
62  "Data in VTK file format are not compressed")
64  "Data in VTK file format are compressed using zlib (not supported yet)")
65  .close();
66 }
67 
68 
69 const int OutputVTK::registrar = Input::register_class< OutputVTK >("vtk") +
71 
72 
73 
75 {
76  this->enable_refinement_ = true;
77 }
78 
79 
80 
82 {
83  this->write_tail();
84 }
85 
86 
87 
88 
90 {
91  ASSERT_PTR(output_mesh_).error();
92 
93  /* It's possible now to do output to the file only in the first process */
94  if(this->rank != 0) {
95  /* TODO: do something, when support for Parallel VTK is added */
96  return 0;
97  }
98 
99  if (! this->_base_file.is_open()) {
100  this->fix_main_file_extension(".pvd");
101  this->_base_filename.open_stream( this->_base_file );
102 
103  LogOut() << "Writing flow output file: " << this->_base_filename << " ... ";
104 
105  this->make_subdirectory();
106  this->write_head();
107  }
108 
109  ostringstream ss;
110  ss << main_output_basename_ << "-"
111  << std::setw(6) << std::setfill('0') << this->current_step
112  << ".vtu";
113 
114 
115  std::string frame_file_name = ss.str();
116  FilePath frame_file_path({main_output_dir_, main_output_basename_, frame_file_name}, FilePath::output_file);
117 
118  /* Set up data file */
119  frame_file_path.open_stream(_data_file);
120 
121  LogOut() << __func__ << ": Writing output file " << this->_base_filename << " ... ";
122 
123  /* Set floating point precision to max */
124  this->_base_file.precision(std::numeric_limits<double>::digits10);
125 
126  /* Strip out relative path and add "base/" string */
127  std::string relative_frame_file = main_output_basename_ + "/" + frame_file_name;
128  this->_base_file << scientific << "<DataSet timestep=\"" << (isfinite(this->time)?this->time:0)
129  << "\" group=\"\" part=\"0\" file=\"" << relative_frame_file <<"\"/>" << endl;
130 
131  LogOut() << "O.K.";
132 
133  LogOut() << __func__ << ": Writing output (frame " << this->current_step << ") file " << relative_frame_file << " ... ";
134 
135  this->write_vtk_vtu();
136 
137  /* Close stream for file of current frame */
138  _data_file.close();
139  //delete data_file;
140  //this->_data_file = NULL;
141 
142  LogOut() << "O.K.";
143 
144  return 1;
145 }
146 
147 
148 
149 
151 {
152  ASSERT_EQ(this->_base_filename.extension(), ".pvd").error();
155 
156  vector<string> sub_path = { main_output_dir_, main_output_basename_, "__tmp__" };
157  FilePath fp(sub_path, FilePath::output_file);
158  fp.create_output_dir();
159 }
160 
161 
162 
163 
165 {
166  ofstream &file = this->_data_file;
167 
168  file << "<?xml version=\"1.0\"?>" << endl;
169  // TODO: test endianess of platform (this would be important, when raw
170  // data will be saved to the VTK file)
171  file << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">" << endl;
172  file << "<UnstructuredGrid>" << endl;
173 }
174 
175 
176 
178 {
179  auto offsets = output_mesh_->offsets_->data_;
180  unsigned int n_elements = offsets.size();
181 
182  data.resize(n_elements);
183  int n_nodes;
184 
185  n_nodes = offsets[0];
186  switch(n_nodes) {
187  case 2:
188  data[0] = (unsigned int)VTK_LINE;
189  break;
190  case 3:
191  data[0] = (unsigned int)VTK_TRIANGLE;
192  break;
193  case 4:
194  data[0] = (unsigned int)VTK_TETRA;
195  break;
196  }
197 
198  for(unsigned int i=1; i < n_elements; i++)
199  {
200  n_nodes = offsets[i]-offsets[i-1];
201  switch(n_nodes) {
202  case 2:
203  data[i] = (unsigned int)VTK_LINE;
204  break;
205  case 3:
206  data[i] = (unsigned int)VTK_TRIANGLE;
207  break;
208  case 4:
209  data[i] = (unsigned int)VTK_TETRA;
210  break;
211  }
212  }
213 }
214 
215 
216 
218 {
219  ofstream &file = this->_data_file;
220 
221  file << "<DataArray type=\"" << vtk_value_type_map(type) << "\" ";
222  // possibly write name
223  if( ! output_data->output_field_name.empty())
224  file << "Name=\"" << output_data->output_field_name <<"\" ";
225  // write number of components
226  if (output_data->n_elem_ > 1)
227  {
228  file
229  << "NumberOfComponents=\"" << output_data->n_elem_ << "\" ";
230  }
231  file << "format=\"ascii\">"
232  << endl;
233 
234  /* Set precision to max */
235  file.precision(std::numeric_limits<double>::digits10);
236 
237  output_data->print_all(file);
238 
239  file << "\n</DataArray>" << endl;
240 
241 }
242 
243 
245 {
246  for(OutputDataPtr data : output_data_vec)
248 }
249 
250 
251 
252 
254  OutputDataFieldVec &output_data_vec)
255 {
256  if (output_data_vec.empty()) return;
257 
258  file << "Scalars=\"";
259  for(OutputDataPtr data : output_data_vec )
260  if (data->n_elem_ == OutputDataBase::N_SCALAR) file << data->output_field_name << ",";
261  file << "\" ";
262 
263  file << "Vectors=\"";
264  for(OutputDataPtr data : output_data_vec )
265  if (data->n_elem_ == OutputDataBase::N_VECTOR) file << data->output_field_name << ",";
266  file << "\" ";
267 
268  file << "Tensors=\"";
269  for(OutputDataPtr data : output_data_vec )
270  if (data->n_elem_ == OutputDataBase::N_TENSOR) file << data->output_field_name << ",";
271  file << "\"";
272 }
273 
274 
276 {
277  ofstream &file = this->_data_file;
278 
279  // merge node and corner data
280  OutputDataFieldVec node_corner_data(output_data_vec_[NODE_DATA]);
281  node_corner_data.insert(node_corner_data.end(),
283 
284  if( ! node_corner_data.empty() ) {
285  /* Write <PointData begin */
286  file << "<PointData ";
287  write_vtk_data_names(file, node_corner_data);
288  file << ">" << endl;
289 
290  /* Write data on nodes */
291  this->write_vtk_data_ascii(output_data_vec_[NODE_DATA]);
292 
293  /* Write data in corners of elements */
295 
296  /* Write PointData end */
297  file << "</PointData>" << endl;
298  }
299 }
300 
301 
303 {
304  ofstream &file = this->_data_file;
305 
306  auto &data_map = this->output_data_vec_[ELEM_DATA];
307  if (data_map.empty()) return;
308 
309  /* Write CellData begin */
310  file << "<CellData ";
311  write_vtk_data_names(file, data_map);
312  file << ">" << endl;
313 
314  /* Write own data */
315  this->write_vtk_data_ascii(data_map);
316 
317  /* Write PointData end */
318  file << "</CellData>" << endl;
319 }
320 
321 
323 {
324  ofstream &file = this->_data_file;
325 
326  file << "</UnstructuredGrid>" << endl;
327  file << "</VTKFile>" << endl;
328 }
329 
330 
332 {
333  ofstream &file = this->_data_file;
334 
335  /* Write header */
336  this->write_vtk_vtu_head();
337 
338  /* When there is no discontinuous data, then write classical vtu */
339  if ( this->output_data_vec_[CORNER_DATA].empty() )
340  {
341  /* Write Piece begin */
342  file << "<Piece NumberOfPoints=\"" << output_mesh_->n_nodes()
343  << "\" NumberOfCells=\"" << output_mesh_->n_elements() <<"\">" << endl;
344 
345  /* Write VTK Geometry */
346  file << "<Points>" << endl;
348  file << "</Points>" << endl;
349 
350 
351  /* Write VTK Topology */
352  file << "<Cells>" << endl;
353  write_vtk_data_ascii(output_mesh_->connectivity_, VTK_INT32 );
355  auto types = std::make_shared<MeshData<unsigned int>>("types");
356  fill_element_types_vector(types->data_);
358  file << "</Cells>" << endl;
359 
360  /* Write VTK scalar and vector data on nodes to the file */
361  this->write_vtk_node_data();
362 
363  /* Write VTK data on elements */
364  this->write_vtk_element_data();
365 
366  /* Write Piece end */
367  file << "</Piece>" << endl;
368 
369  } else {
370  /* Write Piece begin */
371  file << "<Piece NumberOfPoints=\"" << output_mesh_discont_->n_nodes()
372  << "\" NumberOfCells=\"" << output_mesh_->n_elements() <<"\">" << endl;
373 
374  /* Write VTK Geometry */
375  file << "<Points>" << endl;
377  file << "</Points>" << endl;
378 
379  /* Write VTK Topology */
380  file << "<Cells>" << endl;
383  auto types = std::make_shared<MeshData<unsigned int>>("types");
384  fill_element_types_vector(types->data_);
386  file << "</Cells>" << endl;
387 
388  /* Write VTK scalar and vector data on nodes to the file */
389  this->write_vtk_node_data();
390 
391  /* Write VTK data on elements */
392  this->write_vtk_element_data();
393 
394  /* Write Piece end */
395  file << "</Piece>" << endl;
396  }
397 
398  /* Write tail */
399  this->write_vtk_vtu_tail();
400 }
401 
402 
403 
405 {
406  /* It's possible now to do output to the file only in the first process */
407  if(this->rank != 0) {
408  /* TODO: do something, when support for Parallel VTK is added */
409  return 0;
410  }
411 
412  LogOut() << __func__ << ": Writing output file (head) " << this->_base_filename << " ... ";
413 
414  this->_base_file << "<?xml version=\"1.0\"?>" << endl;
415  this->_base_file << "<VTKFile type=\"Collection\" version=\"0.1\" byte_order=\"LittleEndian\">" << endl;
416  this->_base_file << "<Collection>" << endl;
417 
418  LogOut() << "O.K.";
419 
420  return 1;
421 }
422 
423 
425 {
426  /* It's possible now to do output to the file only in the first process */
427  if(this->rank != 0) {
428  /* TODO: do something, when support for Parallel VTK is added */
429  return 0;
430  }
431 
432  LogOut() << __func__ << ": Writing output file (tail) " << this->_base_filename << " ... ";
433 
434  this->_base_file << "</Collection>" << endl;
435  this->_base_file << "</VTKFile>" << endl;
436 
437  LogOut() << "O.K.";
438 
439  return 1;
440 }
441 
442 
443 
444 
Classes for auxiliary output mesh.
string stem() const
Definition: file_path.cc:193
double time
Definition: output_time.hh:208
std::shared_ptr< OutputMesh > output_mesh_
Output mesh.
Definition: output_time.hh:248
std::shared_ptr< OutputMeshDiscontinuous > output_mesh_discont_
Discontinuous (non-conforming) mesh. Used for CORNER_DATA.
Definition: output_time.hh:250
void fix_main_file_extension(std::string extension)
Definition: output_time.cc:157
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:582
void make_subdirectory()
Definition: output_vtk.cc:150
static const Input::Type::Record & get_input_type()
The definition of input record for vtk file format.
Definition: output_vtk.cc:33
void create_output_dir()
Definition: file_path.cc:176
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:331
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:50
Class for declaration of the input of type Bool.
Definition: type_base.hh:434
FilePath _base_filename
Definition: output_time.hh:234
static const int registrar
Registrar of class to factory.
Definition: output_vtk.hh:139
static Input::Type::Abstract & get_input_format_type()
The specification of output file format.
Definition: output_time.cc:58
void fill_element_types_vector(std::vector< unsigned int > &data)
Fills the given vector with VTK element types indicators.
Definition: output_vtk.cc:177
string main_output_basename_
Definition: output_vtk.hh:209
int write_head(void)
This function writes header of VTK (.pvd) file format.
Definition: output_vtk.cc:404
string main_output_dir_
Definition: output_vtk.hh:211
void write_vtk_element_data(void)
Write data on elements to the VTK file (.vtu)
Definition: output_vtk.cc:302
#define LogOut()
Macro defining &#39;log&#39; record of log.
Definition: logger.hh:237
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:286
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:195
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:89
ofstream _data_file
Definition: output_vtk.hh:202
static const Input::Type::Selection & get_input_type_compression()
The definition of input record for selection of compression type.
Definition: output_vtk.cc:59
int current_step
Definition: output_time.hh:203
Selection & add_value(const int value, const std::string &key, const std::string &description="")
Adds one new value with name given by key to the Selection.
void write_vtk_vtu_tail(void)
Write tail of VTK file (.vtu)
Definition: output_vtk.cc:322
This class is used for output data to VTK file format.
Definition: output_vtk.hh:33
string parent_path() const
Definition: file_path.cc:183
static const std::string vtk_value_type_map(VTKValueType t)
Definition: output_vtk.hh:134
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:468
int write_tail(void)
This function writes tail of VTK (.pvd) file format.
Definition: output_vtk.cc:424
ofstream _base_file
Definition: output_time.hh:229
~OutputVTK()
The destructor of this class. It writes tail of the file too.
Definition: output_vtk.cc:81
string extension() const
Definition: file_path.cc:198
Dedicated class for storing path to input and output files.
Definition: file_path.hh:48
void write_vtk_data_ascii(OutputDataFieldVec &output_data_map)
Definition: output_vtk.cc:244
OutputVTK()
The constructor of this class. The head of file is written, when constructor is called.
Definition: output_vtk.cc:74
#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.
void write_vtk_node_data(void)
Write data on nodes to the VTK file (.vtu)
Definition: output_vtk.cc:275
OutputDataFieldVec output_data_vec_[N_DISCRETE_SPACES]
Definition: output_time.hh:198
std::shared_ptr< OutputDataBase > OutputDataPtr
Definition: output_time.hh:191
Record type proxy class.
Definition: type_record.hh:171
bool enable_refinement_
Auxliary flag for refinement enabling, due to gmsh format.
Definition: output_time.hh:255
void write_vtk_vtu_head(void)
Write header of VTK file (.vtu)
Definition: output_vtk.cc:164
Template for classes storing finite set of named values.
#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:49
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:253