Flow123d  build_with_4.0.3-6210fd3
msh_pvdreader.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 msh_pvdreader.cc
15  * @brief
16  * @author dalibor
17  */
18 
19 
20 #include <iostream>
21 #include <vector>
22 #include <pugixml.hpp>
23 
24 #include "msh_pvdreader.hh"
25 #include "msh_vtkreader.hh"
26 
27 
28 
30 : BaseMeshReader(file_name)
31 {
32  data_section_name_ = "DataArray";
33  has_compatible_mesh_ = false;
34  pvd_path_dir_ = file_name.parent_path();
36 }
37 
38 
40 {
41  for (auto file_data : file_list_) {
42  if (file_data.reader != nullptr) delete file_data.reader;
43  }
44 }
45 
46 
48  // will be implemented later
49  // ASSERT_PERMANENT(0).error("Not implemented!");
50 }
51 
52 
54  file_list_[0].reader = new VtkMeshReader(file_list_[0].file_name, this->element_data_values_, file_list_[0].time);
55  file_list_[0].reader->read_nodes(mesh);
56 }
57 
58 
60  file_list_[0].reader->read_elements(mesh);
61 }
62 
63 
65  list_it_->reader->read_element_data(data_cache, header);
66 }
67 
68 
70  pugi::xml_document doc;;
71  doc.load_file( tok_.f_name().c_str() );
72  pugi::xml_node node = doc.child("VTKFile").child("Collection");
73 
74  double last_time = -numeric_limits<double>::infinity(); // check ascending order of times
75  std::vector<std::string> sub_paths; //allow construct paths of VTK files
76  sub_paths.resize(2);
77  sub_paths[0] = pvd_path_dir_;
78 
79  // read PVD data
80  for (pugi::xml_node subnode = node.child("DataSet"); subnode; subnode = subnode.next_sibling("DataSet")) {
81  double time = subnode.attribute("timestep").as_double();
82  if (time <= last_time) {
83  WarningOut().fmt("Wrong time order in PVD file '{}', time '{}'. Skipping this time step.\n", tok_.f_name(), time );
84  } else {
85  sub_paths[1] = subnode.attribute("file").as_string();
86  FilePath vtk_path(sub_paths, FilePath::input_file);
87  last_time = time;
88  file_list_.push_back( VtkFileData(time, vtk_path) );
89  }
90  }
91 }
92 
93 
95  auto comp = [](double t, const VtkFileData &a) {
96  return t * (1.0 + 2.0*numeric_limits<double>::epsilon()) < a.time;
97  };
98 
99  // find iterator to data of VTK file
100  list_it_ = std::upper_bound(file_list_.begin(),
101  file_list_.end(),
102  header_query.time,
103  comp);
104  --list_it_;
105 
106  // check if VTK reader exists and eventually creates its
107  if (!list_it_->reader) {
108  list_it_->reader = new VtkMeshReader(list_it_->file_name, this->element_data_values_, list_it_->time);
109  list_it_->reader->bulk_elements_id_ = this->bulk_elements_id_;
110  list_it_->reader->boundary_elements_id_ = this->boundary_elements_id_;
111  list_it_->reader->has_compatible_mesh_ = true;
112  }
113 
114  return list_it_->reader->find_header(header_query);
115 }
116 
vector< LongIdx > boundary_elements_id_
vector< LongIdx > bulk_elements_id_
Tokenizer tok_
Tokenizer used for reading ASCII file format.
std::string data_section_name_
Store name of field data section specify for type of mesh file.
std::shared_ptr< ElementDataFieldMap > element_data_values_
Cache with last read element data.
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
@ input_file
Definition: file_path.hh:68
string parent_path() const
Definition: file_path.cc:183
Definition: mesh.h:362
std::vector< VtkFileData > file_list_
Store list of VTK files and time steps declared in PVD file.
void read_element_data(ElementDataCacheBase &data_cache, MeshDataHeader header) override
~PvdMeshReader()
Destructor.
void make_header_table() override
std::vector< VtkFileData >::iterator list_it_
Iterator to items of file_list_.
MeshDataHeader & find_header(HeaderQuery &header_query) override
void read_elements(Mesh *mesh) override
void read_physical_names(Mesh *mesh) override
std::string pvd_path_dir_
Path to PVD file allows construct FilePath objects of VTK files.
void read_nodes(Mesh *mesh) override
PvdMeshReader(const FilePath &file_name)
#define WarningOut()
Macro defining 'warning' record of log.
Definition: logger.hh:278
double time
Time of field data (used only for GMSH and PVD reader)
double time
Time of field data (used only for GMSH reader)
Represents data of one VTK file defined in PVD file.