libMesh
unv_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/libmesh_logging.h"
22 #include "libmesh/unv_io.h"
23 #include "libmesh/mesh_base.h"
24 #include "libmesh/edge_edge2.h"
25 #include "libmesh/face_quad4.h"
26 #include "libmesh/face_tri3.h"
27 #include "libmesh/face_tri6.h"
28 #include "libmesh/face_quad8.h"
29 #include "libmesh/face_quad9.h"
30 #include "libmesh/cell_tet4.h"
31 #include "libmesh/cell_hex8.h"
32 #include "libmesh/cell_hex20.h"
33 #include "libmesh/cell_tet10.h"
34 #include "libmesh/cell_prism6.h"
35 #include "libmesh/utility.h"
36 
37 // C++ includes
38 #include <array>
39 #include <iomanip>
40 #include <algorithm> // for std::sort
41 #include <fstream>
42 #include <ctype.h> // isspace
43 #include <sstream> // std::istringstream
44 #include <unordered_map>
45 
46 #ifdef LIBMESH_HAVE_GZSTREAM
47 # include "libmesh/ignore_warnings.h" // shadowing in gzstream.h
48 # include "gzstream.h" // For reading/writing compressed streams
49 # include "libmesh/restore_warnings.h"
50 #endif
51 
52 
53 
54 namespace libMesh
55 {
56 
57 //-----------------------------------------------------------------------------
58 // UNVIO class static members
59 const std::string UNVIO::_nodes_dataset_label = "2411";
60 const std::string UNVIO::_elements_dataset_label = "2412";
61 const std::string UNVIO::_groups_dataset_label = "2467";
62 
63 
64 
65 // ------------------------------------------------------------
66 // UNVIO class members
67 
71  _verbose (false)
72 {
73 }
74 
75 
76 
79  _verbose (false)
80 {
81 }
82 
83 
84 
86 {
87 }
88 
89 
90 
91 bool & UNVIO::verbose ()
92 {
93  return _verbose;
94 }
95 
96 
97 
98 void UNVIO::read (const std::string & file_name)
99 {
100  if (file_name.rfind(".gz") < file_name.size())
101  {
102 #ifdef LIBMESH_HAVE_GZSTREAM
103 
104  igzstream in_stream (file_name.c_str());
105  this->read_implementation (in_stream);
106 
107 #else
108 
109  libmesh_error_msg("ERROR: You must have the zlib.h header files and libraries to read and write compressed streams.");
110 
111 #endif
112  return;
113  }
114 
115  else
116  {
117  std::ifstream in_stream (file_name.c_str());
118  this->read_implementation (in_stream);
119  return;
120  }
121 }
122 
123 
124 void UNVIO::read_implementation (std::istream & in_stream)
125 {
126  // Keep track of what kinds of elements this file contains
127  elems_of_dimension.clear();
128  elems_of_dimension.resize(4, false);
129 
130  {
131  if (!in_stream.good())
132  libmesh_error_msg("ERROR: Input file not good.");
133 
134  // Flags to be set when certain sections are encountered
135  bool
136  found_node = false,
137  found_elem = false,
138  found_group = false;
139 
140  // strings for reading the file line by line
141  std::string
142  old_line,
143  current_line;
144 
145  while (true)
146  {
147  // Save off the old_line. This will provide extra reliability
148  // for detecting the beginnings of the different sections.
149  old_line = current_line;
150 
151  // Try to read something. This may set EOF!
152  std::getline(in_stream, current_line);
153 
154  // If the stream is still "valid", parse the line
155  if (in_stream)
156  {
157  // UNV files always have some amount of leading
158  // whitespace, let's not rely on exactly how much... This
159  // command deletes it.
160  current_line.erase(std::remove_if(current_line.begin(), current_line.end(), isspace),
161  current_line.end());
162 
163  // Parse the nodes section
164  if (current_line == _nodes_dataset_label &&
165  old_line == "-1")
166  {
167  found_node = true;
168  this->nodes_in(in_stream);
169  }
170 
171  // Parse the elements section
172  else if (current_line == _elements_dataset_label &&
173  old_line == "-1")
174  {
175  // The current implementation requires the nodes to
176  // have been read before reaching the elements
177  // section.
178  if (!found_node)
179  libmesh_error_msg("ERROR: The Nodes section must come before the Elements section of the UNV file!");
180 
181  found_elem = true;
182  this->elements_in(in_stream);
183  }
184 
185  // Parse the groups section
186  else if (current_line == _groups_dataset_label &&
187  old_line == "-1")
188  {
189  // The current implementation requires the nodes and
190  // elements to have already been read before reaching
191  // the groups section.
192  if (!found_node || !found_elem)
193  libmesh_error_msg("ERROR: The Nodes and Elements sections must come before the Groups section of the UNV file!");
194 
195  found_group = true;
196  this->groups_in(in_stream);
197  }
198 
199  // We can stop reading once we've found the nodes, elements,
200  // and group sections.
201  if (found_node && found_elem && found_group)
202  break;
203 
204  // If we made it here, we're not done yet, so keep reading
205  continue;
206  }
207 
208  // if (!in_stream) check to see if EOF was set. If so, break out of while loop.
209  if (in_stream.eof())
210  break;
211 
212  // If !in_stream and !in_stream.eof(), stream is in a bad state!
213  libmesh_error_msg("Stream is bad! Perhaps the file does not exist?");
214  } // end while (true)
215 
216  // By now we better have found the datasets for nodes and elements,
217  // otherwise the unv file is bad!
218  if (!found_node)
219  libmesh_error_msg("ERROR: Could not find nodes!");
220 
221  if (!found_elem)
222  libmesh_error_msg("ERROR: Could not find elements!");
223  }
224 
225 
226 
227  {
228  // Set the mesh dimension to the largest element dimension encountered
229  MeshInput<MeshBase>::mesh().set_mesh_dimension(max_elem_dimension_seen());
230 
231 #if LIBMESH_DIM < 3
232  if (MeshInput<MeshBase>::mesh().mesh_dimension() > LIBMESH_DIM)
233  libmesh_error_msg("Cannot open dimension " \
234  << MeshInput<MeshBase>::mesh().mesh_dimension() \
235  << " mesh file when configured without " \
236  << MeshInput<MeshBase>::mesh().mesh_dimension() \
237  << "D support." );
238 #endif
239 
240  // Delete any lower-dimensional elements that might have been
241  // added to the mesh strictly for setting BCs.
242  {
243  // Grab reference to the Mesh
245 
246  unsigned char max_dim = this->max_elem_dimension_seen();
247 
248  for (const auto & elem : mesh.element_ptr_range())
249  if (elem->dim() < max_dim)
250  mesh.delete_elem(elem);
251  }
252 
253  if (this->verbose())
254  libMesh::out << " Finished." << std::endl << std::endl;
255  }
256 }
257 
258 
259 
260 
261 
262 void UNVIO::write (const std::string & file_name)
263 {
264  if (file_name.rfind(".gz") < file_name.size())
265  {
266 #ifdef LIBMESH_HAVE_GZSTREAM
267 
268  ogzstream out_stream(file_name.c_str());
269  this->write_implementation (out_stream);
270 
271 #else
272 
273  libmesh_error_msg("ERROR: You must have the zlib.h header files and libraries to read and write compressed streams.");
274 
275 #endif
276 
277  return;
278  }
279 
280  else
281  {
282  std::ofstream out_stream (file_name.c_str());
283  this->write_implementation (out_stream);
284  return;
285  }
286 }
287 
288 
289 
290 
291 void UNVIO::write_implementation (std::ostream & out_file)
292 {
293  if (!out_file.good())
294  libmesh_error_msg("ERROR: Output file not good.");
295 
296  // write the nodes, then the elements
297  this->nodes_out (out_file);
298  this->elements_out (out_file);
299 }
300 
301 
302 
303 
304 void UNVIO::nodes_in (std::istream & in_file)
305 {
306  LOG_SCOPE("nodes_in()","UNVIO");
307 
308  if (this->verbose())
309  libMesh::out << " Reading nodes" << std::endl;
310 
312 
313  // node label, we use an int here so we can read in a -1
314  int node_label;
315 
316  // We'll just read the floating point values as strings even when
317  // there are no "D" characters in the file. This will make parsing
318  // the file a bit slower but the logic to do so will all be
319  // simplified and in one place...
320  std::string line;
321  std::istringstream coords_stream;
322 
323  // Continue reading nodes until there are none left
324  unsigned ctr = 0;
325  while (true)
326  {
327  // Read the node label
328  in_file >> node_label;
329 
330  // Break out of the while loop when we hit -1
331  if (node_label == -1)
332  break;
333 
334  // Read and discard the the rest of the node data on this line
335  // which we do not currently use:
336  // .) exp_coord_sys_num
337  // .) disp_coord_sys_num
338  // .) color
339  std::getline(in_file, line);
340 
341  // read the floating-point data
342  std::getline(in_file, line);
343 
344  // Replace any "D" characters used for exponents with "E"
345  size_t last_pos = 0;
346  while (true)
347  {
348  last_pos = line.find("D", last_pos);
349 
350  if (last_pos != std::string::npos)
351  line.replace(last_pos, 1, "E");
352  else
353  break;
354  }
355 
356  // Construct a stream object from the current line and then
357  // stream in the xyz values.
358  coords_stream.str(line);
359  coords_stream.clear(); // clear iostate bits! Important!
360 
361  // always 3 coordinates in the UNV file, no matter
362  // what LIBMESH_DIM is.
363  std::array<Real, 3> xyz;
364 
365  coords_stream >> xyz[0] >> xyz[1] >> xyz[2];
366 
367  Point p(xyz[0]);
368 #if LIBMESH_DIM > 1
369  p(1) = xyz[1];
370 #else
371  libmesh_assert_equal_to(xyz[1], 0);
372 #endif
373 #if LIBMESH_DIM > 2
374  p(2) = xyz[2];
375 #else
376  libmesh_assert_equal_to(xyz[2], 0);
377 #endif
378 
379  // Add node to the Mesh, bump the counter.
380  Node * added_node = mesh.add_point(p, ctr++);
381 
382  // Maintain the mapping between UNV node ids and libmesh Node
383  // pointers.
384  _unv_node_id_to_libmesh_node_ptr[node_label] = added_node;
385  }
386 }
387 
388 
389 
391 {
392  unsigned char max_dim = 0;
393 
394  unsigned char elem_dimensions_size = cast_int<unsigned char>
395  (elems_of_dimension.size());
396  // The elems_of_dimension array is 1-based in the UNV reader
397  for (unsigned char i=1; i<elem_dimensions_size; ++i)
398  if (elems_of_dimension[i])
399  max_dim = i;
400 
401  return max_dim;
402 }
403 
404 
405 
406 bool UNVIO::need_D_to_e (std::string & number)
407 {
408  // find "D" in string, start looking at 6th element since dont expect a "D" earlier.
409  std::string::size_type position = number.find("D", 6);
410 
411  if (position != std::string::npos)
412  {
413  // replace "D" in string
414  number.replace(position, 1, "e");
415  return true;
416  }
417  else
418  return false;
419 }
420 
421 
422 
423 void UNVIO::groups_in (std::istream & in_file)
424 {
425  // Grab reference to the Mesh, so we can add boundary info data to it
427 
428  // Record the max and min element dimension seen while reading the file.
429  unsigned char max_dim = this->max_elem_dimension_seen();
430 
431  // Container which stores lower-dimensional elements (based on
432  // Elem::key()) for later assignment of boundary conditions. We
433  // could use e.g. an unordered_set with Elems as keys for this, but
434  // this turns out to be much slower on the search side, since we
435  // have to build an entire side in order to search, rather than just
436  // calling elem->key(side) to compute a key.
437  typedef std::unordered_multimap<dof_id_type, Elem *> map_type;
438  map_type provide_bcs;
439 
440  // Read groups until there aren't any more to read...
441  while (true)
442  {
443  // If we read a -1, it means there is nothing else to read in this section.
444  int group_number;
445  in_file >> group_number;
446 
447  if (group_number == -1)
448  break;
449 
450  // The first record consists of 8 entries:
451  // Field 1 -- group number (that we just read)
452  // Field 2 -- active constraint set no. for group
453  // Field 3 -- active restraint set no. for group
454  // Field 4 -- active load set no. for group
455  // Field 5 -- active dof set no. for group
456  // Field 6 -- active temperature set no. for group
457  // Field 7 -- active contact set no. for group
458  // Field 8 -- number of entities in group
459  // Only the first and last of these are relevant to us...
460  unsigned num_entities;
461  std::string group_name;
462  {
463  unsigned dummy;
464  in_file >> dummy >> dummy >> dummy >> dummy >> dummy >> dummy
465  >> num_entities;
466 
467  // The second record has 1 field, the group name
468  in_file >> group_name;
469  }
470 
471  // The dimension of the elements in the group will determine
472  // whether this is a sideset group or a subdomain group.
473  bool
474  is_subdomain_group = false,
475  is_sideset_group = false;
476 
477  // Each subsequent line has 8 entries, there are two entity type
478  // codes and two tags per line that we need to read. In all my
479  // examples, the entity type code was always "8", which stands for
480  // "finite elements" but it's possible that we could eventually
481  // support other entity type codes...
482  // Field 1 -- entity type code
483  // Field 2 -- entity tag
484  // Field 3 -- entity node leaf id.
485  // Field 4 -- entity component/ ham id.
486  // Field 5 -- entity type code
487  // Field 6 -- entity tag
488  // Field 7 -- entity node leaf id.
489  // Field 8 -- entity component/ ham id.
490  {
491  unsigned entity_type_code, entity_tag, dummy;
492  for (unsigned entity=0; entity<num_entities; ++entity)
493  {
494  in_file >> entity_type_code >> entity_tag >> dummy >> dummy;
495 
496  if (entity_type_code != 8)
497  libMesh::err << "Warning, unrecognized entity type code = "
498  << entity_type_code
499  << " in group "
500  << group_name
501  << std::endl;
502 
503  // Try to find this ID in the map from UNV element ids to libmesh ids.
504  auto it = _unv_elem_id_to_libmesh_elem_id.find(entity_tag);
505  if (it != _unv_elem_id_to_libmesh_elem_id.end())
506  {
507  unsigned libmesh_elem_id = it->second;
508 
509  // Attempt to get a pointer to the elem listed in the group
510  Elem * group_elem = mesh.elem_ptr(libmesh_elem_id);
511 
512  // dim < max_dim means the Elem defines a boundary condition
513  if (group_elem->dim() < max_dim)
514  {
515  is_sideset_group = true;
516 
517  // We can only handle elements that are *exactly*
518  // one dimension lower than the max element
519  // dimension. Not sure if "edge" BCs in 3D
520  // actually make sense/are required...
521  if (group_elem->dim()+1 != max_dim)
522  libmesh_error_msg("ERROR: Expected boundary element of dimension " \
523  << max_dim-1 << " but got " << group_elem->dim());
524 
525  // Set the current group number as the lower-dimensional element's subdomain ID.
526  // We will use this later to set a boundary ID.
527  group_elem->subdomain_id() =
528  cast_int<subdomain_id_type>(group_number);
529 
530  // Store the lower-dimensional element in the provide_bcs container.
531  provide_bcs.insert(std::make_pair(group_elem->key(), group_elem));
532  }
533 
534  // dim == max_dim means this group defines a subdomain ID
535  else if (group_elem->dim() == max_dim)
536  {
537  is_subdomain_group = true;
538  group_elem->subdomain_id() =
539  cast_int<subdomain_id_type>(group_number);
540  }
541 
542  else
543  libmesh_error_msg("ERROR: Found an elem with dim=" \
544  << group_elem->dim() << " > " << "max_dim=" << +max_dim);
545  }
546  else
547  libMesh::err << "WARNING: UNV Element " << entity_tag << " not found while parsing groups." << std::endl;
548  } // end for (entity)
549  } // end scope
550 
551  // Associate this group_number with the group_name in the BoundaryInfo object.
552  if (is_sideset_group)
554  (cast_int<boundary_id_type>(group_number)) = group_name;
555 
556  if (is_subdomain_group)
558  (cast_int<subdomain_id_type>(group_number)) = group_name;
559 
560  } // end while (true)
561 
562  // Loop over elements and try to assign boundary information
563  for (auto & elem : mesh.active_element_ptr_range())
564  if (elem->dim() == max_dim)
565  for (auto sn : elem->side_index_range())
566  for (const auto & pr : as_range(provide_bcs.equal_range (elem->key(sn))))
567  {
568  // This is a max-dimension element that may require BCs.
569  // For each of its sides, including internal sides, we'll
570  // see if a lower-dimensional element provides boundary
571  // information for it. Note that we have not yet called
572  // find_neighbors(), so we can't use elem->neighbor(sn) in
573  // this algorithm...
574 
575  // Build a side to confirm the hash mapped to the correct side.
576  std::unique_ptr<Elem> side (elem->build_side_ptr(sn));
577 
578  // Get a pointer to the lower-dimensional element
579  Elem * lower_dim_elem = pr.second;
580 
581  // This was a hash, so it might not be perfect. Let's verify...
582  if (*lower_dim_elem == *side)
583  mesh.get_boundary_info().add_side(elem, sn, lower_dim_elem->subdomain_id());
584  }
585 }
586 
587 
588 
589 void UNVIO::elements_in (std::istream & in_file)
590 {
591  LOG_SCOPE("elements_in()","UNVIO");
592 
593  if (this->verbose())
594  libMesh::out << " Reading elements" << std::endl;
595 
597 
598  // node label, we use an int here so we can read in a -1
599  int element_label;
600 
601  unsigned
602  n_nodes, // number of nodes on element
603  fe_descriptor_id, // FE descriptor id
604  phys_prop_tab_num, // physical property table number (not supported yet)
605  mat_prop_tab_num, // material property table number (not supported yet)
606  color; // color (not supported yet)
607 
608  // vector that temporarily holds the node labels defining element
609  std::vector<unsigned int> node_labels (21);
610 
611  // vector that assigns element nodes to their correct position
612  // for example:
613  // 44:plane stress | QUAD4
614  // linear quadrilateral |
615  // position in UNV-file | position in libmesh
616  // assign_elem_node[1] = 0
617  // assign_elem_node[2] = 3
618  // assign_elem_node[3] = 2
619  // assign_elem_node[4] = 1
620  //
621  // UNV is 1-based, we leave the 0th element of the vectors unused in order
622  // to prevent confusion, this way we can store elements with up to 20 nodes
623  unsigned int assign_elem_nodes[21];
624 
625  // Read elements until there are none left
626  unsigned ctr = 0;
627  while (true)
628  {
629  // read element label, break out when we read -1
630  in_file >> element_label;
631 
632  if (element_label == -1)
633  break;
634 
635  in_file >> fe_descriptor_id // read FE descriptor id
636  >> phys_prop_tab_num // (not supported yet)
637  >> mat_prop_tab_num // (not supported yet)
638  >> color // (not supported yet)
639  >> n_nodes; // read number of nodes on element
640 
641  // For "beam" type elements, the next three numbers are:
642  // .) beam orientation node number
643  // .) beam fore-end cross section number
644  // .) beam aft-end cross section number
645  // which we discard in libmesh. The "beam" type elements:
646  // 11 Rod
647  // 21 Linear beam
648  // 22 Tapered beam
649  // 23 Curved beam
650  // 24 Parabolic beam
651  // all have fe_descriptor_id < 25.
652  // http://www.sdrl.uc.edu/universal-file-formats-for-modal-analysis-testing-1/file-format-storehouse/unv_2412.htm
653  if (fe_descriptor_id < 25)
654  {
655  unsigned dummy;
656  in_file >> dummy >> dummy >> dummy;
657  }
658 
659  // read node labels (1-based)
660  for (unsigned int j=1; j<=n_nodes; j++)
661  in_file >> node_labels[j];
662 
663  // element pointer, to be allocated
664  Elem * elem = nullptr;
665 
666  switch (fe_descriptor_id)
667  {
668  case 11: // Rod
669  {
670  elem = new Edge2;
671 
672  assign_elem_nodes[1]=0;
673  assign_elem_nodes[2]=1;
674  break;
675  }
676 
677  case 41: // Plane Stress Linear Triangle
678  case 91: // Thin Shell Linear Triangle
679  {
680  elem = new Tri3; // create new element
681 
682  assign_elem_nodes[1]=0;
683  assign_elem_nodes[2]=2;
684  assign_elem_nodes[3]=1;
685  break;
686  }
687 
688  case 42: // Plane Stress Quadratic Triangle
689  case 92: // Thin Shell Quadratic Triangle
690  {
691  elem = new Tri6; // create new element
692 
693  assign_elem_nodes[1]=0;
694  assign_elem_nodes[2]=5;
695  assign_elem_nodes[3]=2;
696  assign_elem_nodes[4]=4;
697  assign_elem_nodes[5]=1;
698  assign_elem_nodes[6]=3;
699  break;
700  }
701 
702  case 43: // Plane Stress Cubic Triangle
703  libmesh_error_msg("ERROR: UNV-element type 43: Plane Stress Cubic Triangle not supported.");
704 
705  case 44: // Plane Stress Linear Quadrilateral
706  case 94: // Thin Shell Linear Quadrilateral
707  {
708  elem = new Quad4; // create new element
709 
710  assign_elem_nodes[1]=0;
711  assign_elem_nodes[2]=3;
712  assign_elem_nodes[3]=2;
713  assign_elem_nodes[4]=1;
714  break;
715  }
716 
717  case 45: // Plane Stress Quadratic Quadrilateral
718  case 95: // Thin Shell Quadratic Quadrilateral
719  {
720  elem = new Quad8; // create new element
721 
722  assign_elem_nodes[1]=0;
723  assign_elem_nodes[2]=7;
724  assign_elem_nodes[3]=3;
725  assign_elem_nodes[4]=6;
726  assign_elem_nodes[5]=2;
727  assign_elem_nodes[6]=5;
728  assign_elem_nodes[7]=1;
729  assign_elem_nodes[8]=4;
730  break;
731  }
732 
733  case 300: // Thin Shell Quadratic Quadrilateral (nine nodes)
734  {
735  elem = new Quad9; // create new element
736 
737  assign_elem_nodes[1]=0;
738  assign_elem_nodes[2]=7;
739  assign_elem_nodes[3]=3;
740  assign_elem_nodes[4]=6;
741  assign_elem_nodes[5]=2;
742  assign_elem_nodes[6]=5;
743  assign_elem_nodes[7]=1;
744  assign_elem_nodes[8]=4;
745  assign_elem_nodes[9]=8;
746  break;
747  }
748 
749  case 46: // Plane Stress Cubic Quadrilateral
750  libmesh_error_msg("ERROR: UNV-element type 46: Plane Stress Cubic Quadrilateral not supported.");
751 
752  case 111: // Solid Linear Tetrahedron
753  {
754  elem = new Tet4; // create new element
755 
756  assign_elem_nodes[1]=0;
757  assign_elem_nodes[2]=1;
758  assign_elem_nodes[3]=2;
759  assign_elem_nodes[4]=3;
760  break;
761  }
762 
763  case 112: // Solid Linear Prism
764  {
765  elem = new Prism6; // create new element
766 
767  assign_elem_nodes[1]=0;
768  assign_elem_nodes[2]=1;
769  assign_elem_nodes[3]=2;
770  assign_elem_nodes[4]=3;
771  assign_elem_nodes[5]=4;
772  assign_elem_nodes[6]=5;
773  break;
774  }
775 
776  case 115: // Solid Linear Brick
777  {
778  elem = new Hex8; // create new element
779 
780  assign_elem_nodes[1]=0;
781  assign_elem_nodes[2]=4;
782  assign_elem_nodes[3]=5;
783  assign_elem_nodes[4]=1;
784  assign_elem_nodes[5]=3;
785  assign_elem_nodes[6]=7;
786  assign_elem_nodes[7]=6;
787  assign_elem_nodes[8]=2;
788  break;
789  }
790 
791  case 116: // Solid Quadratic Brick
792  {
793  elem = new Hex20; // create new element
794 
795  assign_elem_nodes[1]=0;
796  assign_elem_nodes[2]=12;
797  assign_elem_nodes[3]=4;
798  assign_elem_nodes[4]=16;
799  assign_elem_nodes[5]=5;
800  assign_elem_nodes[6]=13;
801  assign_elem_nodes[7]=1;
802  assign_elem_nodes[8]=8;
803 
804  assign_elem_nodes[9]=11;
805  assign_elem_nodes[10]=19;
806  assign_elem_nodes[11]=17;
807  assign_elem_nodes[12]=9;
808 
809  assign_elem_nodes[13]=3;
810  assign_elem_nodes[14]=15;
811  assign_elem_nodes[15]=7;
812  assign_elem_nodes[16]=18;
813  assign_elem_nodes[17]=6;
814  assign_elem_nodes[18]=14;
815  assign_elem_nodes[19]=2;
816  assign_elem_nodes[20]=10;
817  break;
818  }
819 
820  case 117: // Solid Cubic Brick
821  libmesh_error_msg("Error: UNV-element type 117: Solid Cubic Brick not supported.");
822 
823  case 118: // Solid Quadratic Tetrahedron
824  {
825  elem = new Tet10; // create new element
826 
827  assign_elem_nodes[1]=0;
828  assign_elem_nodes[2]=4;
829  assign_elem_nodes[3]=1;
830  assign_elem_nodes[4]=5;
831  assign_elem_nodes[5]=2;
832  assign_elem_nodes[6]=6;
833  assign_elem_nodes[7]=7;
834  assign_elem_nodes[8]=8;
835  assign_elem_nodes[9]=9;
836  assign_elem_nodes[10]=3;
837  break;
838  }
839 
840  default: // Unrecognized element type
841  libmesh_error_msg("ERROR: UNV-element type " << fe_descriptor_id << " not supported.");
842  }
843 
844  // nodes are being stored in element
845  for (dof_id_type j=1; j<=n_nodes; j++)
846  {
847  // Map the UNV node ID to the libmesh node ID
848  auto & node_ptr =
849  libmesh_map_find(_unv_node_id_to_libmesh_node_ptr, node_labels[j]);
850 
851  elem->set_node(assign_elem_nodes[j]) = node_ptr;
852  }
853 
854  elems_of_dimension[elem->dim()] = true;
855 
856  // Set the element's ID
857  elem->set_id(ctr);
858 
859  // Maintain a map from the libmesh (0-based) numbering to the
860  // UNV numbering.
861  _unv_elem_id_to_libmesh_elem_id[element_label] = ctr;
862 
863  // Add the element to the Mesh
864  mesh.add_elem(elem);
865 
866  // Increment the counter for the next iteration
867  ctr++;
868  } // end while (true)
869 }
870 
871 
872 
873 
874 
875 
876 void UNVIO::nodes_out (std::ostream & out_file)
877 {
878  if (this->verbose())
879  libMesh::out << " Writing " << MeshOutput<MeshBase>::mesh().n_nodes() << " nodes" << std::endl;
880 
881  // Write beginning of dataset
882  out_file << " -1\n"
883  << " "
885  << '\n';
886 
887  unsigned int
888  exp_coord_sys_dummy = 0, // export coordinate sys. (not supported yet)
889  disp_coord_sys_dummy = 0, // displacement coordinate sys. (not supp. yet)
890  color_dummy = 0; // color(not supported yet)
891 
892  // A reference to the parent class's mesh
894 
895  // Use scientific notation with capital E and default 17 digits for printing out the coordinates
896  out_file << std::scientific << std::setprecision(this->ascii_precision()) << std::uppercase;
897 
898  for (const auto & current_node : mesh.node_ptr_range())
899  {
900  dof_id_type node_id = current_node->id();
901 
902  out_file << std::setw(10) << node_id
903  << std::setw(10) << exp_coord_sys_dummy
904  << std::setw(10) << disp_coord_sys_dummy
905  << std::setw(10) << color_dummy
906  << '\n';
907 
908  // The coordinates - always write out three coords regardless of LIBMESH_DIM
909  Real x = (*current_node)(0);
910 
911 #if LIBMESH_DIM > 1
912  Real y = (*current_node)(1);
913 #else
914  Real y = 0.;
915 #endif
916 
917 #if LIBMESH_DIM > 2
918  Real z = (*current_node)(2);
919 #else
920  Real z = 0.;
921 #endif
922 
923  out_file << std::setw(25) << x
924  << std::setw(25) << y
925  << std::setw(25) << z
926  << '\n';
927  }
928 
929  // Write end of dataset
930  out_file << " -1\n";
931 }
932 
933 
934 
935 
936 
937 
938 void UNVIO::elements_out(std::ostream & out_file)
939 {
940  if (this->verbose())
941  libMesh::out << " Writing elements" << std::endl;
942 
943  // Write beginning of dataset
944  out_file << " -1\n"
945  << " "
947  << "\n";
948 
949  unsigned
950  fe_descriptor_id = 0, // FE descriptor id
951  phys_prop_tab_dummy = 2, // physical property (not supported yet)
952  mat_prop_tab_dummy = 1, // material property (not supported yet)
953  color_dummy = 7; // color (not supported yet)
954 
955  // vector that assigns element nodes to their correct position
956  // currently only elements with up to 20 nodes
957  //
958  // Example:
959  // QUAD4 | 44:plane stress
960  // | linear quad
961  // position in libMesh | UNV numbering
962  // (note: 0-based) | (note: 1-based)
963  //
964  // assign_elem_node[0] = 0
965  // assign_elem_node[1] = 3
966  // assign_elem_node[2] = 2
967  // assign_elem_node[3] = 1
968  unsigned int assign_elem_nodes[20];
969 
970  unsigned int n_elem_written=0;
971 
972  // A reference to the parent class's mesh
974 
975  for (const auto & elem : mesh.element_ptr_range())
976  {
977  switch (elem->type())
978  {
979 
980  case TRI3:
981  {
982  fe_descriptor_id = 41; // Plane Stress Linear Triangle
983  assign_elem_nodes[0] = 0;
984  assign_elem_nodes[1] = 2;
985  assign_elem_nodes[2] = 1;
986  break;
987  }
988 
989  case TRI6:
990  {
991  fe_descriptor_id = 42; // Plane Stress Quadratic Triangle
992  assign_elem_nodes[0] = 0;
993  assign_elem_nodes[1] = 5;
994  assign_elem_nodes[2] = 2;
995  assign_elem_nodes[3] = 4;
996  assign_elem_nodes[4] = 1;
997  assign_elem_nodes[5] = 3;
998  break;
999  }
1000 
1001  case QUAD4:
1002  {
1003  fe_descriptor_id = 44; // Plane Stress Linear Quadrilateral
1004  assign_elem_nodes[0] = 0;
1005  assign_elem_nodes[1] = 3;
1006  assign_elem_nodes[2] = 2;
1007  assign_elem_nodes[3] = 1;
1008  break;
1009  }
1010 
1011  case QUAD8:
1012  {
1013  fe_descriptor_id = 45; // Plane Stress Quadratic Quadrilateral
1014  assign_elem_nodes[0] = 0;
1015  assign_elem_nodes[1] = 7;
1016  assign_elem_nodes[2] = 3;
1017  assign_elem_nodes[3] = 6;
1018  assign_elem_nodes[4] = 2;
1019  assign_elem_nodes[5] = 5;
1020  assign_elem_nodes[6] = 1;
1021  assign_elem_nodes[7] = 4;
1022  break;
1023  }
1024 
1025  case QUAD9:
1026  {
1027  fe_descriptor_id = 300; // Plane Stress Quadratic Quadrilateral
1028  assign_elem_nodes[0] = 0;
1029  assign_elem_nodes[1] = 7;
1030  assign_elem_nodes[2] = 3;
1031  assign_elem_nodes[3] = 6;
1032  assign_elem_nodes[4] = 2;
1033  assign_elem_nodes[5] = 5;
1034  assign_elem_nodes[6] = 1;
1035  assign_elem_nodes[7] = 4;
1036  assign_elem_nodes[8] = 8;
1037  break;
1038  }
1039 
1040  case TET4:
1041  {
1042  fe_descriptor_id = 111; // Solid Linear Tetrahedron
1043  assign_elem_nodes[0] = 0;
1044  assign_elem_nodes[1] = 1;
1045  assign_elem_nodes[2] = 2;
1046  assign_elem_nodes[3] = 3;
1047  break;
1048  }
1049 
1050  case PRISM6:
1051  {
1052  fe_descriptor_id = 112; // Solid Linear Prism
1053  assign_elem_nodes[0] = 0;
1054  assign_elem_nodes[1] = 1;
1055  assign_elem_nodes[2] = 2;
1056  assign_elem_nodes[3] = 3;
1057  assign_elem_nodes[4] = 4;
1058  assign_elem_nodes[5] = 5;
1059  break;
1060  }
1061 
1062  case HEX8:
1063  {
1064  fe_descriptor_id = 115; // Solid Linear Brick
1065  assign_elem_nodes[0] = 0;
1066  assign_elem_nodes[1] = 4;
1067  assign_elem_nodes[2] = 5;
1068  assign_elem_nodes[3] = 1;
1069  assign_elem_nodes[4] = 3;
1070  assign_elem_nodes[5] = 7;
1071  assign_elem_nodes[6] = 6;
1072  assign_elem_nodes[7] = 2;
1073  break;
1074  }
1075 
1076  case HEX20:
1077  {
1078  fe_descriptor_id = 116; // Solid Quadratic Brick
1079  assign_elem_nodes[ 0] = 0;
1080  assign_elem_nodes[ 1] = 12;
1081  assign_elem_nodes[ 2] = 4;
1082  assign_elem_nodes[ 3] = 16;
1083  assign_elem_nodes[ 4] = 5;
1084  assign_elem_nodes[ 5] = 13;
1085  assign_elem_nodes[ 6] = 1;
1086  assign_elem_nodes[ 7] = 8;
1087 
1088  assign_elem_nodes[ 8] = 11;
1089  assign_elem_nodes[ 9] = 19;
1090  assign_elem_nodes[10] = 17;
1091  assign_elem_nodes[11] = 9;
1092 
1093  assign_elem_nodes[12] = 3;
1094  assign_elem_nodes[13] = 15;
1095  assign_elem_nodes[14] = 7;
1096  assign_elem_nodes[15] = 18;
1097  assign_elem_nodes[16] = 6;
1098  assign_elem_nodes[17] = 14;
1099  assign_elem_nodes[18] = 2;
1100  assign_elem_nodes[19] = 10;
1101 
1102 
1103  break;
1104  }
1105 
1106  case TET10:
1107  {
1108  fe_descriptor_id = 118; // Solid Quadratic Tetrahedron
1109  assign_elem_nodes[0] = 0;
1110  assign_elem_nodes[1] = 4;
1111  assign_elem_nodes[2] = 1;
1112  assign_elem_nodes[3] = 5;
1113  assign_elem_nodes[4] = 2;
1114  assign_elem_nodes[5] = 6;
1115  assign_elem_nodes[6] = 7;
1116  assign_elem_nodes[7] = 8;
1117  assign_elem_nodes[8] = 9;
1118  assign_elem_nodes[9] = 3;
1119  break;
1120  }
1121 
1122  default:
1123  libmesh_error_msg("ERROR: Element type = " << elem->type() << " not supported in " << "UNVIO!");
1124  }
1125 
1126  dof_id_type elem_id = elem->id();
1127 
1128  out_file << std::setw(10) << elem_id // element ID
1129  << std::setw(10) << fe_descriptor_id // type of element
1130  << std::setw(10) << phys_prop_tab_dummy // not supported
1131  << std::setw(10) << mat_prop_tab_dummy // not supported
1132  << std::setw(10) << color_dummy // not supported
1133  << std::setw(10) << elem->n_nodes() // No. of nodes per element
1134  << '\n';
1135 
1136  for (auto j : elem->node_index_range())
1137  {
1138  // assign_elem_nodes[j]-th node: i.e., j loops over the
1139  // libMesh numbering, and assign_elem_nodes[j] over the
1140  // UNV numbering.
1141  const Node * node_in_unv_order = elem->node_ptr(assign_elem_nodes[j]);
1142 
1143  // new record after 8 id entries
1144  if (j==8 || j==16)
1145  out_file << '\n';
1146 
1147  // Write foreign label for this node
1148  dof_id_type node_id = node_in_unv_order->id();
1149 
1150  out_file << std::setw(10) << node_id;
1151  }
1152 
1153  out_file << '\n';
1154 
1155  n_elem_written++;
1156  }
1157 
1158  if (this->verbose())
1159  libMesh::out << " Finished writing " << n_elem_written << " elements" << std::endl;
1160 
1161  // Write end of dataset
1162  out_file << " -1\n";
1163 }
1164 
1165 
1166 
1167 void UNVIO::read_dataset(std::string file_name)
1168 {
1169  std::ifstream in_stream(file_name.c_str());
1170 
1171  if (!in_stream.good())
1172  libmesh_error_msg("Error opening UNV data file.");
1173 
1174  std::string olds, news, dummy;
1175 
1176  while (true)
1177  {
1178  in_stream >> olds >> news;
1179 
1180  // A "-1" followed by a number means the beginning of a dataset.
1181  while (((olds != "-1") || (news == "-1")) && !in_stream.eof())
1182  {
1183  olds = news;
1184  in_stream >> news;
1185  }
1186 
1187  if (in_stream.eof())
1188  break;
1189 
1190  // We only support reading datasets in the "2414" format.
1191  if (news == "2414")
1192  {
1193  // Ignore the rest of the current line and the next two
1194  // lines that contain analysis dataset label and name.
1195  for (unsigned int i=0; i<3; i++)
1196  std::getline(in_stream, dummy);
1197 
1198  // Read the dataset location, where
1199  // 1: Data at nodes
1200  // 2: Data on elements
1201  // Currently only data on nodes is supported.
1202  unsigned int dataset_location;
1203  in_stream >> dataset_location;
1204 
1205  // Currently only nodal datasets are supported.
1206  if (dataset_location != 1)
1207  libmesh_error_msg("ERROR: Currently only Data at nodes is supported.");
1208 
1209  // Ignore the rest of this line and the next five records.
1210  for (unsigned int i=0; i<6; i++)
1211  std::getline(in_stream, dummy);
1212 
1213  // These data are all of no interest to us...
1214  unsigned int
1215  model_type,
1216  analysis_type,
1217  data_characteristic,
1218  result_type;
1219 
1220  // The type of data (complex, real, float, double etc, see
1221  // below).
1222  unsigned int data_type;
1223 
1224  // The number of floating-point values per entity.
1225  unsigned int num_vals_per_node;
1226 
1227  in_stream >> model_type // not used here
1228  >> analysis_type // not used here
1229  >> data_characteristic // not used here
1230  >> result_type // not used here
1231  >> data_type
1232  >> num_vals_per_node;
1233 
1234  // Ignore the rest of record 9, and records 10-13.
1235  for (unsigned int i=0; i<5; i++)
1236  std::getline(in_stream, dummy);
1237 
1238  // Now get the foreign (aka UNV node) node number and
1239  // the respective nodal data.
1240  int f_n_id;
1241  std::vector<Number> values;
1242 
1243  while (true)
1244  {
1245  in_stream >> f_n_id;
1246 
1247  // If -1 then we have reached the end of the dataset.
1248  if (f_n_id == -1)
1249  break;
1250 
1251  // Resize the values vector (usually data in three
1252  // principle directions, i.e. num_vals_per_node = 3).
1253  values.resize(num_vals_per_node);
1254 
1255  // Read the values for the respective node.
1256  for (unsigned int data_cnt=0; data_cnt<num_vals_per_node; data_cnt++)
1257  {
1258  // Check what data type we are reading.
1259  // 2,4: Real
1260  // 5,6: Complex
1261  // other data types are not supported yet.
1262  // As again, these floats may also be written
1263  // using a 'D' instead of an 'e'.
1264  if (data_type == 2 || data_type == 4)
1265  {
1266  std::string buf;
1267  in_stream >> buf;
1268  this->need_D_to_e(buf);
1269 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1270  values[data_cnt] = Complex(std::atof(buf.c_str()), 0.);
1271 #else
1272  values[data_cnt] = std::atof(buf.c_str());
1273 #endif
1274  }
1275 
1276  else if (data_type == 5 || data_type == 6)
1277  {
1278 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1279  Real re_val, im_val;
1280 
1281  std::string buf;
1282  in_stream >> buf;
1283 
1284  if (this->need_D_to_e(buf))
1285  {
1286  re_val = std::atof(buf.c_str());
1287  in_stream >> buf;
1288  this->need_D_to_e(buf);
1289  im_val = std::atof(buf.c_str());
1290  }
1291  else
1292  {
1293  re_val = std::atof(buf.c_str());
1294  in_stream >> im_val;
1295  }
1296 
1297  values[data_cnt] = Complex(re_val,im_val);
1298 #else
1299 
1300  libmesh_error_msg("ERROR: Complex data only supported when libMesh is configured with --enable-complex!");
1301 #endif
1302  }
1303 
1304  else
1305  libmesh_error_msg("ERROR: Data type not supported.");
1306 
1307  } // end loop data_cnt
1308 
1309  // Get a pointer to the Node associated with the UNV node id.
1310  auto & node_ptr = libmesh_map_find(_unv_node_id_to_libmesh_node_ptr, f_n_id);
1311 
1312  // Store the nodal values in our map which uses the
1313  // libMesh Node* as the key. We use operator[] here
1314  // because we want to create an empty vector if the
1315  // entry does not already exist.
1316  _node_data[node_ptr] = values;
1317  } // end while (true)
1318  } // end if (news == "2414")
1319  } // end while (true)
1320 }
1321 
1322 
1323 
1324 const std::vector<Number> *
1325 UNVIO::get_data (Node * node) const
1326 {
1327  auto it = _node_data.find(node);
1328 
1329  if (it == _node_data.end())
1330  return nullptr;
1331  else
1332  return &(it->second);
1333 }
1334 
1335 
1336 } // namespace libMesh
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::UNVIO::write
virtual void write(const std::string &) override
This method implements writing a mesh to a specified file.
Definition: unv_io.C:262
libMesh::UNVIO::max_elem_dimension_seen
unsigned char max_elem_dimension_seen()
Definition: unv_io.C:390
libMesh::Quad9
The QUAD9 is an element in 2D composed of 9 nodes.
Definition: face_quad9.h:51
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::Tri6
The Tri6 is an element in 2D composed of 6 nodes.
Definition: face_tri6.h:56
libMesh::MeshBase::delete_elem
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
libMesh::MeshBase::active_element_ptr_range
virtual SimpleRange< element_iterator > active_element_ptr_range()=0
libMesh::DofObject::set_id
dof_id_type & set_id()
Definition: dof_object.h:776
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::TET10
Definition: enum_elem_type.h:46
libMesh::BoundaryInfo::sideset_name
std::string & sideset_name(boundary_id_type id)
Definition: boundary_info.C:2341
libMesh::UNVIO::read
virtual void read(const std::string &) override
This method implements reading a mesh from a specified file.
Definition: unv_io.C:98
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::Hex20
The Hex20 is an element in 3D composed of 20 nodes.
Definition: cell_hex20.h:68
libMesh::Prism6
The Prism6 is an element in 3D composed of 6 nodes.
Definition: cell_prism6.h:52
libMesh::UNVIO::elements_in
void elements_in(std::istream &in_file)
Method reads elements and stores them in std::vector<Elem *> _elements in the same order as they come...
Definition: unv_io.C:589
libMesh::MeshBase::elem_ptr
virtual const Elem * elem_ptr(const dof_id_type i) const =0
libMesh::UNVIO::_unv_elem_id_to_libmesh_elem_id
std::map< unsigned, unsigned > _unv_elem_id_to_libmesh_elem_id
Map UNV element IDs to libmesh element IDs.
Definition: unv_io.h:208
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::MeshOutput< MeshBase >::ascii_precision
unsigned int & ascii_precision()
Return/set the precision to use when writing ASCII files.
Definition: mesh_output.h:257
libMesh::UNVIO::_unv_node_id_to_libmesh_node_ptr
std::map< dof_id_type, Node * > _unv_node_id_to_libmesh_node_ptr
Maps UNV node IDs to libMesh Node*s.
Definition: unv_io.h:188
libMesh::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
libMesh::UNVIO::_elements_dataset_label
static const std::string _elements_dataset_label
label for the element dataset
Definition: unv_io.h:198
libMesh::UNVIO::~UNVIO
virtual ~UNVIO()
Destructor.
Definition: unv_io.C:85
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::UNVIO::verbose
bool & verbose()
Set the flag indicating if we should be verbose.
Definition: unv_io.C:91
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::UNVIO::_nodes_dataset_label
static const std::string _nodes_dataset_label
label for the node dataset
Definition: unv_io.h:193
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::UNVIO::read_dataset
void read_dataset(std::string file_name)
Read a UNV data file containing a dataset of type "2414".
Definition: unv_io.C:1167
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::UNVIO::nodes_out
void nodes_out(std::ostream &out_file)
Outputs nodes to the file out_file.
Definition: unv_io.C:876
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::Tet4
The Tet4 is an element in 3D composed of 4 nodes.
Definition: cell_tet4.h:53
libMesh::as_range
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
libMesh::UNVIO::get_data
const std::vector< Number > * get_data(Node *node) const
Definition: unv_io.C:1325
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::Quad4
The QUAD4 is an element in 2D composed of 4 nodes.
Definition: face_quad4.h:51
libMesh::UNVIO::elements_out
void elements_out(std::ostream &out_file)
Outputs the element data to the file out_file.
Definition: unv_io.C:938
libMesh::UNVIO::groups_in
void groups_in(std::istream &in_file)
Reads the "groups" section of the file.
Definition: unv_io.C:423
libMesh::UNVIO::_groups_dataset_label
static const std::string _groups_dataset_label
label for the groups dataset
Definition: unv_io.h:203
libMesh::Elem::set_node
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2059
libMesh::Complex
std::complex< Real > Complex
Definition: libmesh_common.h:160
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::MeshOutput
This class defines an abstract interface for Mesh output.
Definition: mesh_output.h:53
libMesh::UNVIO::write_implementation
void write_implementation(std::ostream &out_stream)
The actual implementation of the write function.
Definition: unv_io.C:291
libMesh::UNVIO::_verbose
bool _verbose
should be be verbose?
Definition: unv_io.h:182
libMesh::Hex8
The Hex8 is an element in 3D composed of 8 nodes.
Definition: cell_hex8.h:53
libMesh::MeshOutput::mesh
const MT & mesh() const
Definition: mesh_output.h:247
libMesh::Tet10
The Tet10 is an element in 3D composed of 10 nodes.
Definition: cell_tet10.h:60
libMesh::Quad8
The QUAD8 is an element in 2D composed of 8 nodes.
Definition: face_quad8.h:51
libMesh::MeshBase::add_elem
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
libMesh::Elem::subdomain_id
subdomain_id_type subdomain_id() const
Definition: elem.h:2069
libMesh::DofObject::id
dof_id_type id() const
Definition: dof_object.h:767
libMesh::UNVIO::UNVIO
UNVIO(MeshBase &mesh)
Constructor.
Definition: unv_io.C:68
libMesh::Edge2
The Edge2 is an element in 1D composed of 2 nodes.
Definition: edge_edge2.h:43
libMesh::Tri3
The Tri3 is an element in 2D composed of 3 nodes.
Definition: face_tri3.h:56
libMesh::MeshBase::subdomain_name
std::string & subdomain_name(subdomain_id_type id)
Definition: mesh_base.C:717
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::UNVIO::read_implementation
void read_implementation(std::istream &in_stream)
The actual implementation of the read function.
Definition: unv_io.C:124
libMesh::MeshInput::mesh
MT & mesh()
Definition: mesh_input.h:169
libMesh::QUAD9
Definition: enum_elem_type.h:43
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::err
OStreamProxy err
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Elem::key
virtual dof_id_type key(const unsigned int s) const =0
libMesh::UNVIO::_node_data
std::map< Node *, std::vector< Number > > _node_data
Map from libMesh Node* to data at that node, as read in by the read_dataset() function.
Definition: unv_io.h:214
libMesh::out
OStreamProxy out
libMesh::UNVIO::nodes_in
void nodes_in(std::istream &in_file)
Read nodes from file.
Definition: unv_io.C:304
libMesh::BoundaryInfo::add_side
void add_side(const dof_id_type elem, const unsigned short int side, const boundary_id_type id)
Add side side of element number elem with boundary id id to the boundary information data structure.
Definition: boundary_info.C:886
libMesh::MeshInput
This class defines an abstract interface for Mesh input.
Definition: mesh_base.h:60
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::QUAD8
Definition: enum_elem_type.h:42
libMesh::UNVIO::need_D_to_e
bool need_D_to_e(std::string &number)
Replaces "1.1111D+00" with "1.1111e+00" if necessary.
Definition: unv_io.C:406