Flow123d  jenkins-Flow123d-linux-release-multijob-282
output_vtk.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
5  *
6  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
7  * especially for academic research:
8  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
9  *
10  * This program is free software; you can redistribute it and/or modify it under the terms
11  * of the GNU General Public License version 3 as published by the Free Software Foundation.
12  *
13  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
14  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
15  * See the GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License along with this program; if not,
18  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
19  *
20  *
21  * $Id: output_vtk.cc 2505 2013-09-13 14:52:27Z jiri.hnidek $
22  * $Revision$
23  * $LastChangedBy$
24  * $LastChangedDate$
25  *
26  * @file output_vtk.cc
27  * @brief The functions for outputs to VTK files.
28  *
29  */
30 
31 #include "output_vtk.hh"
32 #include <limits.h>
33 #include "mesh/mesh.h"
34 #include "output_data_base.hh"
35 
36 
37 using namespace Input::Type;
38 
40  = Record("vtk", "Parameters of vtk output format.")
41  // It is derived from abstract class
43  .declare_key("variant", input_type_variant, Default("ascii"),
44  "Variant of output stream file format.")
45  // The parallel or serial variant
46  .declare_key("parallel", Bool(), Default("false"),
47  "Parallel or serial version of file format.")
48  // Type of compression
49  .declare_key("compression", input_type_compression, Default("none"),
50  "Compression used in output stream file format.");
51 
52 
54  = Selection("VTK variant (ascii or binary)")
56  "ASCII variant of VTK file format")
58  "Binary variant of VTK file format (not supported yet)");
59 
60 
62  = Selection("Type of compression of VTK file format")
64  "Data in VTK file format are not compressed")
66  "Data in VTK file format are compressed using zlib (not supported yet)");
67 
68 
70 {
71  this->fix_main_file_extension(".pvd");
72 
73  if(this->rank == 0) {
74  this->_base_file.open(this->_base_filename.c_str());
75  INPUT_CHECK( this->_base_file.is_open() , "Can not open output file: %s\n", this->_base_filename.c_str() );
76  xprintf(MsgLog, "Writing flow output file: %s ... \n", this->_base_filename.c_str());
77  }
78 
79  this->make_subdirectory();
80  this->write_head();
81 }
82 
83 
84 
86 {}
87 
88 
89 
91 {
92  this->write_tail();
93 }
94 
95 
96 
97 
99 {
100  ASSERT(_mesh != nullptr, "Null mesh.\n");
101 
102  /* It's possible now to do output to the file only in the first process */
103  if(this->rank != 0) {
104  /* TODO: do something, when support for Parallel VTK is added */
105  return 0;
106  }
107 
108  ostringstream ss;
109  ss << main_output_basename_ << "-"
110  << std::setw(6) << std::setfill('0') << this->current_step
111  << ".vtu";
112 
113 
114  std::string frame_file_name = ss.str();
115  std::string frame_file_path = main_output_dir_ + "/" + main_output_basename_ + "/" + frame_file_name;
116 
117  _data_file.open(frame_file_path);
118  if(_data_file.is_open() == false) {
119  xprintf(Err, "Could not write output to the file: %s\n", frame_file_path.c_str());
120  return 0;
121  } else {
122  /* Set up data file */
123 
124  xprintf(MsgLog, "%s: Writing output file %s ... ",
125  __func__, this->_base_filename.c_str());
126 
127 
128  /* Set floating point precision to max */
129  this->_base_file.precision(std::numeric_limits<double>::digits10);
130 
131  /* Strip out relative path and add "base/" string */
132  std::string relative_frame_file = main_output_basename_ + "/" + frame_file_name;
133  this->_base_file << scientific << "<DataSet timestep=\"" << (isfinite(this->time)?this->time:0)
134  << "\" group=\"\" part=\"0\" file=\"" << relative_frame_file <<"\"/>" << endl;
135 
136  xprintf(MsgLog, "O.K.\n");
137 
138  xprintf(MsgLog, "%s: Writing output (frame %d) file %s ... ", __func__,
139  this->current_step, relative_frame_file.c_str());
140 
141  this->write_vtk_vtu();
142 
143  /* Close stream for file of current frame */
144  _data_file.close();
145  //delete data_file;
146  //this->_data_file = NULL;
147 
148  xprintf(MsgLog, "O.K.\n");
149  }
150 
151  return 1;
152 }
153 
154 
155 
156 
158 {
159  string main_file="./" + this->_base_filename; // guarantee that find_last_of succeeds
160  ASSERT( main_file.substr( main_file.size() - 4) == ".pvd" , "none");
161  unsigned int last_sep_pos=main_file.find_last_of(DIR_DELIMITER);
162  main_output_dir_=main_file.substr(2, last_sep_pos-2);
163  main_output_basename_=main_file.substr(last_sep_pos+1);
164  main_output_basename_=main_output_basename_.substr(0, main_output_basename_.size() - 4); // 5 = ".pvd".size() +1
165 
167  fp.create_output_dir();
168 }
169 
170 
171 
172 
173 
175 {
176  ofstream &file = this->_data_file;
177 
178  file << "<?xml version=\"1.0\"?>" << endl;
179  // TODO: test endianess of platform (this would be important, when raw
180  // data will be saved to the VTK file)
181  file << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">" << endl;
182  file << "<UnstructuredGrid>" << endl;
183 }
184 
186 {
187  Mesh *mesh = this->_mesh;
188  ofstream &file = this->_data_file;
189 
190  int tmp;
191 
192  /* Write Points begin*/
193  file << "<Points>" << endl;
194  /* Write DataArray begin */
195  file << "<DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"ascii\">" << endl;
196  /* Write own coordinates */
197  tmp = 0;
198  /* Set floating point precision */
199  file.precision(std::numeric_limits<double>::digits10);
200  FOR_NODES(mesh, node) {
201  node->aux = tmp; /* store index in the auxiliary variable */
202 
203  file << scientific << node->getX() << " ";
204  file << scientific << node->getY() << " ";
205  file << scientific << node->getZ() << " ";
206 
207  tmp++;
208  }
209  /* Write DataArray end */
210  file << endl << "</DataArray>" << endl;
211  /* Write Points end*/
212  file << "</Points>" << endl;
213 }
214 
216 {
217  Mesh *mesh = this->_mesh;
218  ofstream &file = this->_data_file;
219 
220  Node* node;
221  //ElementIter ele;
222  unsigned int li;
223  int tmp;
224 
225  /* Write Cells begin*/
226  file << "<Cells>" << endl;
227  /* Write DataArray begin */
228  file << "<DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">" << endl;
229  /* Write own coordinates */
230  FOR_ELEMENTS(mesh, ele) {
231  FOR_ELEMENT_NODES(ele, li) {
232  node = ele->node[li];
233  file << node->aux << " "; /* Write connectivity */
234  }
235  }
236  /* Write DataArray end */
237  file << endl << "</DataArray>" << endl;
238 
239  /* Write DataArray begin */
240  file << "<DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">" << endl;
241  /* Write number of nodes for each element */
242  tmp = 0;
243  FOR_ELEMENTS(mesh, ele) {
244  switch(ele->dim()) {
245  case 1:
246  tmp += VTK_LINE_SIZE;
247  break;
248  case 2:
249  tmp += VTK_TRIANGLE_SIZE;
250  break;
251  case 3:
252  tmp += VTK_TETRA_SIZE;
253  break;
254  }
255  file << tmp << " ";
256  }
257  /* Write DataArray end */
258  file << endl << "</DataArray>" << endl;
259 
260  /* Write DataArray begin */
261  file << "<DataArray type=\"UInt8\" Name=\"types\" format=\"ascii\">" << endl;
262  /* Write type of nodes for each element */
263  FOR_ELEMENTS(mesh, ele) {
264  switch(ele->dim()) {
265  case 1:
266  file << (int)VTK_LINE << " ";
267  break;
268  case 2:
269  file << (int)VTK_TRIANGLE << " ";
270  break;
271  case 3:
272  file << (int)VTK_TETRA << " ";
273  break;
274  }
275  }
276  /* Write DataArray end */
277  file << endl << "</DataArray>" << endl;
278 
279  /* Write Cells end*/
280  file << "</Cells>" << endl;
281 }
282 
284 {
285  Mesh *mesh = this->_mesh;
286  ofstream &file = this->_data_file;
287 
288  NodeIter node;
289  unsigned int li;
290 
291  /* Write Points begin*/
292  file << "<Points>" << endl;
293  /* Write DataArray begin */
294  file << "<DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"ascii\">" << endl;
295  /* Set floating point precision */
296  file.precision(std::numeric_limits<double>::digits10);
297  FOR_ELEMENTS(mesh, ele) {
298  FOR_ELEMENT_NODES(ele, li) {
299  node = ele->node[li];
300 
301  file << scientific << node->getX() << " ";
302  file << scientific << node->getY() << " ";
303  file << scientific << node->getZ() << " ";
304  }
305  }
306  /* Write DataArray end */
307  file << endl << "</DataArray>" << endl;
308  /* Write Points end*/
309  file << "</Points>" << endl;
310 }
311 
313 {
314  Mesh *mesh = this->_mesh;
315  ofstream &file = this->_data_file;
316 
317  //Node* node;
318  //ElementIter ele;
319  unsigned int li, tmp;
320 
321  /* Write Cells begin*/
322  file << "<Cells>" << endl;
323  /* Write DataArray begin */
324  file << "<DataArray type=\"Int32\" Name=\"connectivity\" format=\"ascii\">" << endl;
325  /* Write own coordinates */
326  tmp = 0;
327  FOR_ELEMENTS(mesh, ele) {
328  FOR_ELEMENT_NODES(ele, li) {
329  file << tmp << " "; /* Write connectivity */
330  tmp++;
331  }
332  }
333  /* Write DataArray end */
334  file << endl << "</DataArray>" << endl;
335 
336  /* Write DataArray begin */
337  file << "<DataArray type=\"Int32\" Name=\"offsets\" format=\"ascii\">" << endl;
338  /* Write number of nodes for each element */
339  tmp = 0;
340  FOR_ELEMENTS(mesh, ele) {
341  switch(ele->dim()) {
342  case 1:
343  tmp += VTK_LINE_SIZE;
344  break;
345  case 2:
346  tmp += VTK_TRIANGLE_SIZE;
347  break;
348  case 3:
349  tmp += VTK_TETRA_SIZE;
350  break;
351  }
352  file << tmp << " ";
353  }
354  /* Write DataArray end */
355  file << endl << "</DataArray>" << endl;
356 
357  /* Write DataArray begin */
358  file << "<DataArray type=\"UInt8\" Name=\"types\" format=\"ascii\">" << endl;
359  /* Write type of nodes for each element */
360  FOR_ELEMENTS(mesh, ele) {
361  switch(ele->dim()) {
362  case 1:
363  file << (int)VTK_LINE << " ";
364  break;
365  case 2:
366  file << (int)VTK_TRIANGLE << " ";
367  break;
368  case 3:
369  file << (int)VTK_TETRA << " ";
370  break;
371  }
372  }
373  /* Write DataArray end */
374  file << endl << "</DataArray>" << endl;
375 
376  /* Write Cells end*/
377  file << "</Cells>" << endl;
378 }
379 
380 
381 
382 
384 {
385  ofstream &file = this->_data_file;
386 
387  for(OutputDataPtr data : output_data_vec)
388  {
389  file << "<DataArray type=\"Float64\" "
390  << "Name=\"" << data->output_field_name <<"\" ";
391  if (data->n_elem_ > 1)
392  {
393  file
394  << "NumberOfComponents=\"" << data->n_elem_ << "\" ";
395  }
396  file << "format=\"ascii\">"
397  << endl;
398 
399  /* Set precision to max */
400  file.precision(std::numeric_limits<double>::digits10);
401 
402  data->print_all(file);
403 
404  file << "\n</DataArray>" << endl;
405  }
406 }
407 
408 
409 
410 
412  OutputDataFieldVec &output_data_vec)
413 {
414  if (output_data_vec.empty()) return;
415 
416  file << "Scalars=\"";
417  for(OutputDataPtr data : output_data_vec )
418  if (data->n_elem_ == OutputDataBase::N_SCALAR) file << data->output_field_name << ",";
419  file << "\" ";
420 
421  file << "Vectors=\"";
422  for(OutputDataPtr data : output_data_vec )
423  if (data->n_elem_ == OutputDataBase::N_VECTOR) file << data->output_field_name << ",";
424  file << "\" ";
425 
426  file << "Tensors=\"";
427  for(OutputDataPtr data : output_data_vec )
428  if (data->n_elem_ == OutputDataBase::N_TENSOR) file << data->output_field_name << ",";
429  file << "\"";
430 }
431 
432 
434 {
435  ofstream &file = this->_data_file;
436 
437  // merge node and corner data
438  OutputDataFieldVec node_corner_data(output_data_vec_[NODE_DATA]);
439  node_corner_data.insert(node_corner_data.end(),
441 
442  if( ! node_corner_data.empty() ) {
443  /* Write <PointData begin */
444  file << "<PointData ";
445  write_vtk_data_names(file, node_corner_data);
446  file << ">" << endl;
447 
448  /* Write data on nodes */
449  this->write_vtk_data_ascii(output_data_vec_[NODE_DATA]);
450 
451  /* Write data in corners of elements */
453 
454  /* Write PointData end */
455  file << "</PointData>" << endl;
456  }
457 }
458 
459 
461 {
462  ofstream &file = this->_data_file;
463 
464  auto &data_map = this->output_data_vec_[ELEM_DATA];
465  if (data_map.empty()) return;
466 
467  /* Write CellData begin */
468  file << "<CellData ";
469  write_vtk_data_names(file, data_map);
470  file << ">" << endl;
471 
472  /* Write own data */
473  this->write_vtk_data_ascii(data_map);
474 
475  /* Write PointData end */
476  file << "</CellData>" << endl;
477 }
478 
479 
481 {
482  ofstream &file = this->_data_file;
483 
484  file << "</UnstructuredGrid>" << endl;
485  file << "</VTKFile>" << endl;
486 }
487 
488 
490 {
491  ofstream &file = this->_data_file;
492  Mesh *mesh = this->_mesh;
493 
494  /* Write header */
495  this->write_vtk_vtu_head();
496 
497  /* When there is no discontinuous data, then write classical vtu */
498  if ( this->output_data_vec_[CORNER_DATA].empty() )
499  {
500  /* Write Piece begin */
501  file << "<Piece NumberOfPoints=\"" << mesh->n_nodes() << "\" NumberOfCells=\"" << mesh->n_elements() <<"\">" << endl;
502 
503  /* Write VTK Geometry */
504  this->write_vtk_geometry();
505 
506  /* Write VTK Topology */
507  this->write_vtk_topology();
508 
509  /* Write VTK scalar and vector data on nodes to the file */
510  this->write_vtk_node_data();
511 
512  /* Write VTK data on elements */
513  this->write_vtk_element_data();
514 
515  /* Write Piece end */
516  file << "</Piece>" << endl;
517 
518  } else {
519  /* Write Piece begin */
520  file << "<Piece NumberOfPoints=\"" << mesh->n_corners() << "\" NumberOfCells=\"" << mesh->n_elements() <<"\">" << endl;
521 
522  /* Write VTK Geometry */
524 
525  /* Write VTK Topology */
527 
528  /* Write VTK scalar and vector data on nodes to the file */
529  this->write_vtk_node_data();
530 
531  /* Write VTK data on elements */
532  this->write_vtk_element_data();
533 
534  /* Write Piece end */
535  file << "</Piece>" << endl;
536  }
537 
538  /* Write tail */
539  this->write_vtk_vtu_tail();
540 }
541 
542 
543 
545 {
546  /* It's possible now to do output to the file only in the first process */
547  if(this->rank != 0) {
548  /* TODO: do something, when support for Parallel VTK is added */
549  return 0;
550  }
551 
552  xprintf(MsgLog, "%s: Writing output file (head) %s ... ", __func__,
553  this->_base_filename.c_str() );
554 
555  this->_base_file << "<?xml version=\"1.0\"?>" << endl;
556  this->_base_file << "<VTKFile type=\"Collection\" version=\"0.1\" byte_order=\"LittleEndian\">" << endl;
557  this->_base_file << "<Collection>" << endl;
558 
559  xprintf(MsgLog, "O.K.\n");
560 
561  return 1;
562 }
563 
564 
566 {
567  /* It's possible now to do output to the file only in the first process */
568  if(this->rank != 0) {
569  /* TODO: do something, when support for Parallel VTK is added */
570  return 0;
571  }
572 
573  xprintf(MsgLog, "%s: Writing output file (tail) %s ... ", __func__,
574  this->_base_filename.c_str() );
575 
576  this->_base_file << "</Collection>" << endl;
577  this->_base_file << "</VTKFile>" << endl;
578 
579  xprintf(MsgLog, "O.K.\n");
580 
581  return 1;
582 }
583 
584 
585 
586 
587 
double time
Definition: output_time.hh:190
Mesh * _mesh
Definition: output_time.hh:222
static Input::Type::AbstractRecord input_format_type
The specification of output file format.
Definition: output_time.hh:62
void write_vtk_discont_topology(void)
Write topology (connection of nodes) to the VTK file (.vtu)
Definition: output_vtk.cc:312
#define FOR_ELEMENT_NODES(i, j)
Definition: elements.h:150
void fix_main_file_extension(std::string extension)
Definition: output_time.cc:110
double getY() const
Definition: nodes.hh:71
void make_subdirectory()
Definition: output_vtk.cc:157
static Input::Type::Selection input_type_variant
The definition of input record for selection of variant of file format.
Definition: output_vtk.hh:70
Definition: nodes.hh:44
void create_output_dir()
Definition: file_path.cc:128
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:489
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:39
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:364
Class for declaration of the input of type Bool.
Definition: type_base.hh:332
???
Definition: mesh.h:109
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:197
void write_vtk_geometry(void)
Write geometry (position of nodes) to the VTK file (.vtu)
Definition: output_vtk.cc:185
string main_output_basename_
Definition: output_vtk.hh:221
int write_head(void)
This function writes header of VTK (.pvd) file format.
Definition: output_vtk.cc:544
string main_output_dir_
Definition: output_vtk.hh:223
void write_vtk_element_data(void)
Write data on elements to the VTK file (.vtu)
Definition: output_vtk.cc:460
Definition: system.hh:72
void write_vtk_discont_geometry(void)
Write geometry (position of nodes) to the VTK file (.vtu)
Definition: output_vtk.cc:283
static Input::Type::Record input_type
The definition of input record for vtk file format.
Definition: output_vtk.hh:65
double getZ() const
Definition: nodes.hh:73
int write_data(void)
This function write data to VTK (.pvd) file format for curent time.
Definition: output_vtk.cc:98
unsigned int n_elements() const
Definition: mesh.h:141
ofstream _data_file
Definition: output_vtk.hh:214
#define ASSERT(...)
Definition: global_defs.h:121
int current_step
Definition: output_time.hh:185
static Input::Type::Selection input_type_compression
The definition of input record for selection of compression type.
Definition: output_vtk.hh:75
Selection & add_value(const int value, const std::string &key, const std::string &description="")
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
#define xprintf(...)
Definition: system.hh:100
double getX() const
Definition: nodes.hh:69
void write_vtk_vtu_tail(void)
Write tail of VTK file (.vtu)
Definition: output_vtk.cc:480
unsigned int n_corners()
Definition: mesh.cc:181
int write_tail(void)
This function writes tail of VTK (.pvd) file format.
Definition: output_vtk.cc:565
int aux
Definition: nodes.hh:105
ofstream _base_file
Definition: output_time.hh:212
The class for outputting data during time.
Definition: output_time.hh:32
~OutputVTK()
The destructor of this class. It writes tail of the file too.
Definition: output_vtk.cc:90
#define FOR_NODES(_mesh_, i)
Definition: mesh.h:72
string _base_filename
Definition: output_time.hh:217
#define INPUT_CHECK(i,...)
Debugging macros.
Definition: global_defs.h:61
Dedicated class for storing path to input and output files.
Definition: file_path.hh:32
void write_vtk_data_ascii(OutputDataFieldVec &output_data_map)
Definition: output_vtk.cc:383
OutputVTK()
The constructor of this class. The head of file is written, when constructor is called.
Definition: output_vtk.cc:85
void write_vtk_node_data(void)
Write data on nodes to the VTK file (.vtu)
Definition: output_vtk.cc:433
#define DIR_DELIMITER
Definition: system.hh:56
OutputDataFieldVec output_data_vec_[N_DISCRETE_SPACES]
Definition: output_time.hh:180
Definition: system.hh:72
unsigned int n_nodes() const
Definition: mesh.h:137
std::shared_ptr< OutputDataBase > OutputDataPtr
Definition: output_time.hh:173
Record type proxy class.
Definition: type_record.hh:169
void write_vtk_vtu_head(void)
Write header of VTK file (.vtu)
Definition: output_vtk.cc:174
Template for classes storing finite set of named values.
void write_vtk_topology(void)
Write topology (connection of nodes) to the VTK file (.vtu)
Definition: output_vtk.cc:215
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:411
Record & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:430