Flow123d  release_2.2.0-20-gb8056ca
msh_vtkreader.hh
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_vtkreader.hh
15  * @brief
16  * @author dalibor
17  */
18 
19 #ifndef MSH_VTK_READER_HH
20 #define MSH_VTK_READER_HH
21 
22 #include <string>
23 #include <istream>
24 #include <pugixml.hpp>
25 
26 #include "io/msh_basereader.hh"
27 #include "system/file_path.hh"
28 
29 class VtkMeshReader : public BaseMeshReader {
30 public:
31  TYPEDEF_ERR_INFO(EI_VTKFile, std::string);
32  TYPEDEF_ERR_INFO(EI_ExpectedFormat, std::string);
33  TYPEDEF_ERR_INFO(EI_ErrMessage, std::string);
34  TYPEDEF_ERR_INFO(EI_SectionTypeName, std::string);
35  TYPEDEF_ERR_INFO(EI_TagType, std::string);
36  TYPEDEF_ERR_INFO(EI_TagName, std::string);
37  DECLARE_EXCEPTION(ExcInvalidFormat,
38  << "Invalid format of DataArray " << EI_FieldName::val << ", expected " << EI_ExpectedFormat::val << "\n"
39  << "in the input file: " << EI_VTKFile::qval);
40  DECLARE_EXCEPTION(ExcUnknownFormat,
41  << "Unsupported or missing format of DataArray " << EI_FieldName::val << "\n" << "in the input file: " << EI_VTKFile::qval);
42  DECLARE_EXCEPTION(ExcWrongType,
43  << EI_ErrMessage::val << " data type of " << EI_SectionTypeName::val << "\n" << "in the input file: " << EI_VTKFile::qval);
44  DECLARE_EXCEPTION(ExcIncompatibleMesh,
45  << "Incompatible meshes, " << EI_ErrMessage::val << "\n" << "for VTK input file: " << EI_VTKFile::qval);
46  DECLARE_EXCEPTION(ExcMissingTag,
47  << "Missing " << EI_TagType::val << " " << EI_TagName::val << "\n" << " in the input file: " << EI_VTKFile::qval);
48 
49  /// Possible data sections in UnstructuredGrid - Piece node.
50  enum DataSections {
52  };
53 
54  /// Type of data formats - ascii, binary or compressed with zLib.
55  enum DataFormat {
57  };
58 
59  /**
60  * Map of DataArray sections in VTK file.
61  *
62  * For each field_name contains MeshDataHeader.
63  */
65 
66  /**
67  * Construct the VTK format reader from given FilePath.
68  * This opens the file for reading.
69  */
70  VtkMeshReader(const FilePath &file_name);
71 
72  /// Destructor
74 
75  /**
76  * Read regions from the VTK file and save the physical sections as regions in the RegionDB.
77  *
78  * Region Labels starting with '!' are treated as boundary regions. Elements of these regions are used just to
79  * assign regions to the boundary and are not used in actual FEM computations.
80  */
81  void read_physical_names(Mesh * mesh) override;
82 
83  /**
84  * Check if nodes and elements of VTK mesh is compatible with \p mesh.
85  *
86  * - to all nodes of VTK mesh must exists one and only one nodes in second mesh
87  * - the same must occur for elements
88  * - method fill vector \p bulk_elements_id_
89  * - it is necessary to call this method before calling \p get_element_data
90  */
91  void check_compatible_mesh(Mesh &mesh) override;
92 
93 protected:
94  /**
95  * private method for reading of nodes
96  */
97  void read_nodes(Mesh * mesh);
98 
99  /**
100  * Method for reading of elements.
101  * Input of the mesh allows changing regions within the input file.
102  */
103  void read_elements(Mesh * mesh);
104 
105  /**
106  * Find header of DataArray section of VTK file given by field_name.
107  *
108  * Note: \p time has no effect (it is only for continuity with GMSH reader).
109  */
110  MeshDataHeader & find_header(double time, std::string field_name) override;
111 
112  /// Reads table of DataArray headers through pugixml interface
113  void make_header_table() override;
114 
115  /// Helper method that create DataArray header of given xml node (used from \p make_header_table)
116  MeshDataHeader create_header(pugi::xml_node node, unsigned int n_entities, Tokenizer::Position pos);
117 
118  /// Get DataType by value of string
119  DataType get_data_type(std::string type_str);
120 
121  /// Return size of value of data_type.
122  unsigned int type_value_size(DataType data_type);
123 
124  /// Parse ascii data to data cache
125  void parse_ascii_data(ElementDataCacheBase &data_cache, unsigned int n_components, unsigned int n_entities,
126  Tokenizer::Position pos, bool boundary_domain);
127 
128  /// Parse binary data to data cache
129  void parse_binary_data(ElementDataCacheBase &data_cache, unsigned int n_components, unsigned int n_entities,
130  Tokenizer::Position pos, bool boundary_domain, DataType value_type);
131 
132  /// Uncompress and parse binary compressed data to data cache
133  void parse_compressed_data(ElementDataCacheBase &data_cache, unsigned int n_components, unsigned int n_entities,
134  Tokenizer::Position pos, bool boundary_domain, DataType value_type);
135 
136  /// Set base attributes of VTK and get count of nodes and elements.
137  void read_base_vtk_attributes(pugi::xml_node vtk_node, unsigned int &n_nodes, unsigned int &n_elements);
138 
139  /// Get position of AppendedData tag in VTK file
140  Tokenizer::Position get_appended_position();
141 
142  /**
143  * Implements @p BaseMeshReader::read_element_data.
144  */
145  void read_element_data(ElementDataCacheBase &data_cache, MeshDataHeader actual_header, unsigned int n_components,
146  bool boundary_domain) override;
147 
148  /**
149  * Compare two points representing by armadillo vector.
150  *
151  * - used in \p check_compatible_mesh method
152  * - calculate with \p point_tolerance parameter
153  */
154  bool compare_points(arma::vec3 &p1, arma::vec3 &p2);
155 
156  /// Tolerance during comparison point data with GMSH nodes.
157  static const double point_tolerance;
158 
159  /// header type of VTK file (only for appended data)
161 
162  /// variants of data format (ascii, appended, compressed appended)
164 
165  /// Table with data of DataArray headers
166  HeaderTable header_table_;
167 
168  /// input stream allow read appended data, used only if this tag exists
169  std::istream *data_stream_;
170 
171  /// store count of read entities
172  unsigned int n_read_;
173 
174 };
175 
176 #endif /* MSH_VTK_READER_HH */
177 
void make_header_table() override
Reads table of DataArray headers through pugixml interface.
bool compare_points(arma::vec3 &p1, arma::vec3 &p2)
DataSections
Possible data sections in UnstructuredGrid - Piece node.
HeaderTable header_table_
Table with data of DataArray headers.
void read_nodes(Mesh *mesh)
static const double point_tolerance
Tolerance during comparison point data with GMSH nodes.
DataFormat data_format_
variants of data format (ascii, appended, compressed appended)
MeshDataHeader & find_header(double time, std::string field_name) override
VtkMeshReader(const FilePath &file_name)
void read_base_vtk_attributes(pugi::xml_node vtk_node, unsigned int &n_nodes, unsigned int &n_elements)
Set base attributes of VTK and get count of nodes and elements.
Definition: mesh.h:97
DataType get_data_type(std::string type_str)
Get DataType by value of string.
void parse_compressed_data(ElementDataCacheBase &data_cache, unsigned int n_components, unsigned int n_entities, Tokenizer::Position pos, bool boundary_domain, DataType value_type)
Uncompress and parse binary compressed data to data cache.
void read_physical_names(Mesh *mesh) override
std::istream * data_stream_
input stream allow read appended data, used only if this tag exists
void check_compatible_mesh(Mesh &mesh) override
TYPEDEF_ERR_INFO(EI_VTKFile, std::string)
~VtkMeshReader()
Destructor.
DECLARE_EXCEPTION(ExcInvalidFormat,<< "Invalid format of DataArray "<< EI_FieldName::val<< ", expected "<< EI_ExpectedFormat::val<< "\n"<< "in the input file: "<< EI_VTKFile::qval)
DataType header_type_
header type of VTK file (only for appended data)
unsigned int n_read_
store count of read entities
Dedicated class for storing path to input and output files.
Definition: file_path.hh:54
void read_elements(Mesh *mesh)
void parse_binary_data(ElementDataCacheBase &data_cache, unsigned int n_components, unsigned int n_entities, Tokenizer::Position pos, bool boundary_domain, DataType value_type)
Parse binary data to data cache.
MeshDataHeader create_header(pugi::xml_node node, unsigned int n_entities, Tokenizer::Position pos)
Helper method that create DataArray header of given xml node (used from make_header_table) ...
std::map< std::string, MeshDataHeader > HeaderTable
void read_element_data(ElementDataCacheBase &data_cache, MeshDataHeader actual_header, unsigned int n_components, bool boundary_domain) override
DataFormat
Type of data formats - ascii, binary or compressed with zLib.
DataType
Types of VTK data (value &#39;undefined&#39; for empty value)
unsigned int type_value_size(DataType data_type)
Return size of value of data_type.
Tokenizer::Position get_appended_position()
Get position of AppendedData tag in VTK file.
void parse_ascii_data(ElementDataCacheBase &data_cache, unsigned int n_components, unsigned int n_entities, Tokenizer::Position pos, bool boundary_domain)
Parse ascii data to data cache.