libMesh
ucd_io.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // Local includes
20 #include "libmesh/libmesh_config.h"
21 #include "libmesh/ucd_io.h"
22 #include "libmesh/mesh_base.h"
23 #include "libmesh/face_quad4.h"
24 #include "libmesh/face_tri3.h"
25 #include "libmesh/cell_tet4.h"
26 #include "libmesh/cell_hex8.h"
27 #include "libmesh/cell_prism6.h"
28 #include "libmesh/enum_io_package.h"
29 #include "libmesh/enum_elem_type.h"
30 #include "libmesh/int_range.h"
31 #include "libmesh/utility.h"
32 
33 #ifdef LIBMESH_HAVE_GZSTREAM
34 # include "libmesh/ignore_warnings.h" // shadowing in gzstream.h
35 # include "gzstream.h" // For reading/writing compressed streams
36 # include "libmesh/restore_warnings.h"
37 #endif
38 
39 // C++ includes
40 #include <array>
41 #include <fstream>
42 
43 
44 namespace libMesh
45 {
46 
47 // Initialize the static data members by calling the static build functions.
48 std::map<ElemType, std::string> UCDIO::_writing_element_map = UCDIO::build_writing_element_map();
49 std::map<std::string, ElemType> UCDIO::_reading_element_map = UCDIO::build_reading_element_map();
50 
51 
52 
53 // Static function used to build the _writing_element_map.
54 std::map<ElemType, std::string> UCDIO::build_writing_element_map()
55 {
56  std::map<ElemType, std::string> ret;
57  ret[EDGE2] = "edge";
58  ret[TRI3] = "tri";
59  ret[QUAD4] = "quad";
60  ret[TET4] = "tet";
61  ret[HEX8] = "hex";
62  ret[PRISM6] = "prism";
63  ret[PYRAMID5] = "pyramid";
64  return ret;
65 }
66 
67 
68 
69 // Static function used to build the _reading_element_map.
70 std::map<std::string, ElemType> UCDIO::build_reading_element_map()
71 {
72  std::map<std::string, ElemType> ret;
73  ret["edge"] = EDGE2;
74  ret["tri"] = TRI3;
75  ret["quad"] = QUAD4;
76  ret["tet"] = TET4;
77  ret["hex"] = HEX8;
78  ret["prism"] = PRISM6;
79  ret["pyramid"] = PYRAMID5;
80  return ret;
81 }
82 
83 
84 void UCDIO::read (const std::string & file_name)
85 {
86  if (file_name.rfind(".gz") < file_name.size())
87  {
88 #ifdef LIBMESH_HAVE_GZSTREAM
89  igzstream in_stream (file_name.c_str());
90  this->read_implementation (in_stream);
91 #else
92  libmesh_error_msg("ERROR: You must have the zlib.h header files and libraries to read and write compressed streams.");
93 #endif
94  }
95 
96  else
97  {
98  std::ifstream in_stream (file_name.c_str());
99  this->read_implementation (in_stream);
100  }
101 }
102 
103 
104 
105 void UCDIO::write (const std::string & file_name)
106 {
107  if (file_name.rfind(".gz") < file_name.size())
108  {
109 #ifdef LIBMESH_HAVE_GZSTREAM
110  ogzstream out_stream (file_name.c_str());
111  this->write_implementation (out_stream);
112 #else
113  libmesh_error_msg("ERROR: You must have the zlib.h header files and libraries to read and write compressed streams.");
114 #endif
115  }
116 
117  else
118  {
119  std::ofstream out_stream (file_name.c_str());
120  this->write_implementation (out_stream);
121  }
122 }
123 
124 
125 
126 void UCDIO::read_implementation (std::istream & in)
127 {
128  // This is a serial-only process for now;
129  // the Mesh should be read on processor 0 and
130  // broadcast later
131  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().processor_id(), 0);
132 
133  // Check input buffer
134  libmesh_assert (in.good());
135 
137 
138  // Keep track of what kinds of elements this file contains
139  elems_of_dimension.clear();
140  elems_of_dimension.resize(4, false);
141 
142  this->skip_comment_lines (in, '#');
143 
144  unsigned int nNodes=0, nElem=0, dummy=0;
145 
146  in >> nNodes // Read the number of nodes from the stream
147  >> nElem // Read the number of elements from the stream
148  >> dummy
149  >> dummy
150  >> dummy;
151 
152 
153  // Read the nodal coordinates. Note that UCD format always
154  // stores (x,y,z), and in 2D z=0. We don't need to store this,
155  // however. So, we read in x,y,z for each node and make a point
156  // in the proper way based on what dimension we're in
157  {
158  for (unsigned int i=0; i<nNodes; i++)
159  {
160  libmesh_assert (in.good());
161 
162  std::array<Real, 3> xyz;
163 
164  in >> dummy // Point number
165  >> xyz[0] // x-coordinate value
166  >> xyz[1] // y-coordinate value
167  >> xyz[2]; // z-coordinate value
168 
169  Point p(xyz[0]);
170 #if LIBMESH_DIM > 1
171  p(1) = xyz[1];
172 #else
173  libmesh_assert_equal_to(xyz[1], 0);
174 #endif
175 #if LIBMESH_DIM > 2
176  p(2) = xyz[2];
177 #else
178  libmesh_assert_equal_to(xyz[2], 0);
179 #endif
180 
181  // Build the node
182  mesh.add_point (p, i);
183  }
184  }
185 
186  // Read the elements from the stream. Notice that the UCD node-numbering
187  // scheme is 1-based, and we just created a 0-based scheme above
188  // (which is of course what we want). So, when we read in the nodal
189  // connectivity for each element we need to take 1 off the value of
190  // each node so that we get the right thing.
191  {
192  unsigned int material_id=0, node=0;
193  std::string type;
194 
195  for (unsigned int i=0; i<nElem; i++)
196  {
197  libmesh_assert (in.good());
198 
199  // The cell type can be either tri, quad, tet, hex, or prism.
200  in >> dummy // Cell number, means nothing to us
201  >> material_id // We'll use this for the element subdomain id.
202  >> type; // string describing cell type
203 
204  // Convert the UCD type string to a libmesh ElemType
205  ElemType et = libmesh_map_find(_reading_element_map, type);
206 
207  // Build the required type and release it into our custody.
208  Elem * elem = Elem::build(et).release();
209 
210  for (auto n : elem->node_index_range())
211  {
212  libmesh_assert (in.good());
213 
214  in >> node; // read the current node
215  node -= 1; // UCD is 1-based, so subtract
216 
217  libmesh_assert_less (node, mesh.n_nodes());
218 
219  // assign the node
220  elem->set_node(n) = mesh.node_ptr(node);
221  }
222 
223  elems_of_dimension[elem->dim()] = true;
224 
225  // Set the element's subdomain ID based on the material_id.
226  elem->subdomain_id() = cast_int<subdomain_id_type>(material_id);
227 
228  // Add the element to the mesh
229  elem->set_id(i);
230  mesh.add_elem (elem);
231  }
232 
233  // Set the mesh dimension to the largest encountered for an element
234  for (unsigned char i=0; i!=4; ++i)
235  if (elems_of_dimension[i])
237 
238 #if LIBMESH_DIM < 3
239  if (mesh.mesh_dimension() > LIBMESH_DIM)
240  libmesh_error_msg("Cannot open dimension " \
241  << mesh.mesh_dimension() \
242  << " mesh file when configured without " \
243  << mesh.mesh_dimension() \
244  << "D support.");
245 #endif
246  }
247 }
248 
249 
250 
251 void UCDIO::write_implementation (std::ostream & out_stream)
252 {
253  libmesh_assert (out_stream.good());
254 
256 
257  // UCD doesn't work any dimension except 3?
258  if (mesh.mesh_dimension() != 3)
259  libmesh_error_msg("Error: Can't write boundary elements for meshes of dimension less than 3. " \
260  << "Mesh dimension = " << mesh.mesh_dimension());
261 
262  // Write header
263  this->write_header(out_stream, mesh, mesh.n_elem(), 0);
264 
265  // Write the node coordinates
266  this->write_nodes(out_stream, mesh);
267 
268  // Write the elements
269  this->write_interior_elems(out_stream, mesh);
270 }
271 
272 
273 
274 void UCDIO::write_header(std::ostream & out_stream,
275  const MeshBase & mesh,
276  dof_id_type n_elems,
277  unsigned int n_vars)
278 {
279  libmesh_assert (out_stream.good());
280  // TODO: We could print out the libmesh revision used to write this file here.
281  out_stream << "# For a description of the UCD format see the AVS Developer's guide.\n"
282  << "#\n";
283 
284  // Write the mesh info
285  out_stream << mesh.n_nodes() << " "
286  << n_elems << " "
287  << n_vars << " "
288  << " 0 0\n";
289 }
290 
291 
292 
293 void UCDIO::write_nodes(std::ostream & out_stream,
294  const MeshBase & mesh)
295 {
296  // 1-based node number for UCD
297  unsigned int n=1;
298 
299  // Write the node coordinates
300  for (auto & node : mesh.node_ptr_range())
301  {
302  libmesh_assert (out_stream.good());
303 
304  out_stream << n++ << "\t";
305  node->write_unformatted(out_stream);
306  }
307 }
308 
309 
310 
311 void UCDIO::write_interior_elems(std::ostream & out_stream,
312  const MeshBase & mesh)
313 {
314  // 1-based element number for UCD
315  unsigned int e=1;
316 
317  // Write element information
318  for (const auto & elem : mesh.element_ptr_range())
319  {
320  libmesh_assert (out_stream.good());
321 
322  // Look up the corresponding UCD element type in the static map.
323  std::string elem_string = libmesh_map_find(_writing_element_map, elem->type());
324 
325  // Write the element's subdomain ID as the UCD "material_id".
326  out_stream << e++ << " " << elem->subdomain_id() << " " << elem_string << "\t";
327  elem->write_connectivity(out_stream, UCD);
328  }
329 }
330 
331 
332 
333 void UCDIO::write_nodal_data(const std::string & fname,
334  const std::vector<Number> & soln,
335  const std::vector<std::string> & names)
336 {
338 
339  const dof_id_type n_elem = mesh.n_elem();
340 
341  // Only processor 0 does the writing
342  if (mesh.processor_id())
343  return;
344 
345  std::ofstream out_stream(fname.c_str());
346 
347  // UCD doesn't work in 1D
349 
350  // Write header
351  this->write_header(out_stream,mesh,n_elem,
352  cast_int<unsigned int>(names.size()));
353 
354  // Write the node coordinates
355  this->write_nodes(out_stream, mesh);
356 
357  // Write the elements
358  this->write_interior_elems(out_stream, mesh);
359 
360  // Write the solution
361  this->write_soln(out_stream, mesh, names, soln);
362 }
363 
364 
365 
366 void UCDIO::write_soln(std::ostream & out_stream,
367  const MeshBase & mesh,
368  const std::vector<std::string> & names,
369  const std::vector<Number> & soln)
370 {
371  libmesh_assert (out_stream.good());
372 
373  // First write out how many variables and how many components per variable
374  out_stream << names.size();
375  for (std::size_t i = 0, ns = names.size(); i < ns; i++)
376  {
377  libmesh_assert (out_stream.good());
378  // Each named variable has only 1 component
379  out_stream << " 1";
380  }
381  out_stream << std::endl;
382 
383  // Now write out variable names and units. Since we don't store units
384  // We just write out dummy.
385  for (const auto & name : names)
386  {
387  libmesh_assert (out_stream.good());
388  out_stream << name << ", dummy" << std::endl;
389  }
390 
391  // Now, for each node, write out the solution variables.
392  // We use a 1-based node numbering for UCD.
393  std::size_t nv = names.size();
394  for (auto n : IntRange<std::size_t>(1, mesh.n_nodes()))
395  {
396  libmesh_assert (out_stream.good());
397  out_stream << n;
398 
399  for (std::size_t var = 0; var != nv; var++)
400  {
401  std::size_t idx = nv*(n-1) + var;
402 
403  out_stream << " " << soln[idx];
404  }
405  out_stream << std::endl;
406  }
407 }
408 
409 } // namespace libMesh
libMesh::UCDIO::_writing_element_map
static std::map< ElemType, std::string > _writing_element_map
Definition: ucd_io.h:149
libMesh::UCDIO::write_implementation
void write_implementation(std::ostream &out_stream)
The actual implementation of the write function.
Definition: ucd_io.C:251
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::MeshTools::n_elem
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Count up the number of elements of a specific type (as defined by an iterator range).
Definition: mesh_tools.C:705
n_vars
unsigned int n_vars
Definition: adaptivity_ex3.C:116
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::MeshInput< MeshBase >::skip_comment_lines
void skip_comment_lines(std::istream &in, const char comment_start)
Reads input from in, skipping all the lines that start with the character comment_start.
Definition: mesh_input.h:179
libMesh::MeshBase::n_elem
virtual dof_id_type n_elem() const =0
libMesh::DofObject::set_id
dof_id_type & set_id()
Definition: dof_object.h:776
libMesh::UCDIO::read_implementation
void read_implementation(std::istream &in_stream)
The actual implementation of the read function.
Definition: ucd_io.C:126
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Elem::dim
virtual unsigned short dim() const =0
libMesh::Elem::node_index_range
IntRange< unsigned short > node_index_range() const
Definition: elem.h:2170
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::node_ptr
virtual const Node * node_ptr(const dof_id_type i) const =0
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::UCDIO::write_soln
void write_soln(std::ostream &out, const MeshBase &mesh, const std::vector< std::string > &names, const std::vector< Number > &soln)
Writes all nodal solution variables.
Definition: ucd_io.C:366
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::MeshBase::node_ptr_range
virtual SimpleRange< node_iterator > node_ptr_range()=0
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::MeshTools::Generation::Private::idx
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
A useful inline function which replaces the macros used previously.
Definition: mesh_generation.C:72
libMesh::UCDIO::build_reading_element_map
static std::map< std::string, ElemType > build_reading_element_map()
Definition: ucd_io.C:70
libMesh::UCDIO::_reading_element_map
static std::map< std::string, ElemType > _reading_element_map
Definition: ucd_io.h:153
libMesh::Elem::set_node
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2059
libMesh::UCDIO::read
virtual void read(const std::string &) override
This method implements reading a mesh from a specified file in UCD format.
Definition: ucd_io.C:84
libMesh::MeshBase::n_nodes
virtual dof_id_type n_nodes() const =0
libMesh::MeshOutput
This class defines an abstract interface for Mesh output.
Definition: mesh_output.h:53
libMesh::MeshOutput::mesh
const MT & mesh() const
Definition: mesh_output.h:247
libMesh::MeshBase::add_elem
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
libMesh::PYRAMID5
Definition: enum_elem_type.h:53
libMesh::UCDIO::write_header
void write_header(std::ostream &out, const MeshBase &mesh, dof_id_type n_elems, unsigned int n_vars)
Write UCD format header.
Definition: ucd_io.C:274
libMesh::Elem::subdomain_id
subdomain_id_type subdomain_id() const
Definition: elem.h:2069
libMesh::UCDIO::write
virtual void write(const std::string &) override
This method implements writing a mesh to a specified file in UCD format.
Definition: ucd_io.C:105
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::UCDIO::build_writing_element_map
static std::map< ElemType, std::string > build_writing_element_map()
Definition: ucd_io.C:54
libMesh::UCDIO::write_nodal_data
virtual void write_nodal_data(const std::string &fname, const std::vector< Number > &soln, const std::vector< std::string > &names) override
This method implements writing a mesh and solution to a specified file in UCD format.
Definition: ucd_io.C:333
libMesh::MeshInput< MeshBase >::mesh
MeshBase & mesh()
Definition: mesh_input.h:169
libMesh::MeshBase::add_point
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
libMesh::UCDIO::write_interior_elems
void write_interior_elems(std::ostream &out, const MeshBase &mesh)
Write element information.
Definition: ucd_io.C:311
libMesh::MeshBase::set_mesh_dimension
void set_mesh_dimension(unsigned char d)
Resets the logical dimension of the mesh.
Definition: mesh_base.h:218
libMesh::Elem::build
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:246
libMesh::UCD
Definition: enum_io_package.h:45
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
libMesh::UCDIO::write_nodes
void write_nodes(std::ostream &out, const MeshBase &mesh)
Write node information.
Definition: ucd_io.C:293
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::MeshInput< MeshBase >::elems_of_dimension
std::vector< bool > elems_of_dimension
A vector of bools describing what dimension elements have been encountered when reading a mesh.
Definition: mesh_input.h:97
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33