libMesh
gmsh_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 // C++ includes
19 #include <fstream>
20 #include <set>
21 #include <cstring> // std::memcpy
22 #include <numeric>
23 #include <unordered_map>
24 #include <cstddef>
25 
26 // Local includes
27 #include "libmesh/libmesh_config.h"
28 #include "libmesh/libmesh_logging.h"
29 #include "libmesh/boundary_info.h"
30 #include "libmesh/elem.h"
31 #include "libmesh/gmsh_io.h"
32 #include "libmesh/mesh_base.h"
33 #include "libmesh/int_range.h"
34 #include "libmesh/utility.h" // map_find
35 
36 namespace libMesh
37 {
38 
39 // Initialize the static data member
41 
42 
43 
44 // Definition of the static function which constructs the ElementMaps object.
46 {
47  // Object to be filled up
48  ElementMaps em;
49 
50  // POINT (import only)
51  em.in.insert(std::make_pair(15, ElementDefinition(NODEELEM, 15, 0, 1)));
52 
53  // Add elements with trivial node mappings
54  em.add_def(ElementDefinition(EDGE2, 1, 1, 2));
55  em.add_def(ElementDefinition(EDGE3, 8, 1, 3));
56  em.add_def(ElementDefinition(TRI3, 2, 2, 3));
57  em.add_def(ElementDefinition(TRI6, 9, 2, 6));
58  em.add_def(ElementDefinition(QUAD4, 3, 2, 4));
59  em.add_def(ElementDefinition(QUAD8, 16, 2, 8));
60  em.add_def(ElementDefinition(QUAD9, 10, 2, 9));
61  em.add_def(ElementDefinition(HEX8, 5, 3, 8));
62  em.add_def(ElementDefinition(TET4, 4, 3, 4));
63  em.add_def(ElementDefinition(PRISM6, 6, 3, 6));
64  em.add_def(ElementDefinition(PYRAMID5, 7, 3, 5));
65 
66  // Add elements with non-trivial node mappings
67 
68  // HEX20
69  {
70  ElementDefinition eledef(HEX20, 17, 3, 20);
71  const unsigned int nodes[] = {0,1,2,3,4,5,6,7,8,11,12,9,13,10,14,15,16,19,17,18};
72  std::vector<unsigned int>(nodes, nodes+eledef.nnodes).swap(eledef.nodes); // swap trick
73  em.add_def(eledef);
74  }
75 
76  // HEX27
77  {
78  ElementDefinition eledef(HEX27, 12, 3, 27);
79  const unsigned int nodes[] = {0,1,2,3,4,5,6,7,8,11,12,9,13,10,14,
80  15,16,19,17,18,20,21,24,22,23,25,26};
81  std::vector<unsigned int>(nodes, nodes+eledef.nnodes).swap(eledef.nodes); // swap trick
82  em.add_def(eledef);
83  }
84 
85  // TET10
86  {
87  ElementDefinition eledef(TET10, 11, 3, 10);
88  const unsigned int nodes[] = {0,1,2,3,4,5,6,7,9,8};
89  std::vector<unsigned int>(nodes, nodes+eledef.nnodes).swap(eledef.nodes); // swap trick
90  em.add_def(eledef);
91  }
92 
93  // PRISM15
94  {
95  ElementDefinition eledef(PRISM15, 18, 3, 15);
96  const unsigned int nodes[] = {0,1,2,3,4,5,6,8,9,7,10,11,12,14,13};
97  std::vector<unsigned int>(nodes, nodes+eledef.nnodes).swap(eledef.nodes); // swap trick
98  em.add_def(eledef);
99  }
100 
101  // PRISM18
102  {
103  ElementDefinition eledef(PRISM18, 13, 3, 18);
104  const unsigned int nodes[] = {0,1,2,3,4,5,6,8,9,7,10,11,12,14,13,15,17,16};
105  std::vector<unsigned int>(nodes, nodes+eledef.nnodes).swap(eledef.nodes); // swap trick
106  em.add_def(eledef);
107  }
108 
109  return em;
110 }
111 
112 
113 
116  _binary(false),
117  _write_lower_dimensional_elements(true)
118 {
119 }
120 
121 
122 
126  _binary (false),
127  _write_lower_dimensional_elements(true)
128 {
129 }
130 
131 
132 
133 bool & GmshIO::binary ()
134 {
135  return _binary;
136 }
137 
138 
139 
141 {
143 }
144 
145 
146 
147 void GmshIO::read (const std::string & name)
148 {
149  std::ifstream in (name.c_str());
150  this->read_mesh (in);
151 }
152 
153 
154 
155 void GmshIO::read_mesh(std::istream & in)
156 {
157  // This is a serial-only process for now;
158  // the Mesh should be read on processor 0 and
159  // broadcast later
160  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().processor_id(), 0);
161 
162  libmesh_assert(in.good());
163 
164  // clear any data in the mesh
166  mesh.clear();
167 
168  // some variables
169  int format=0, size=0;
170  Real version = 1.0;
171 
172  // Keep track of lower-dimensional blocks which are not BCs, but
173  // actually blocks of lower-dimensional elements.
174  std::set<subdomain_id_type> lower_dimensional_blocks;
175 
176  // Mapping from physical id -> (physical dim, physical name) pairs.
177  // These can refer to either "sidesets" or "subdomains"; we need to
178  // wait until the Mesh has been read to know which is which. Note
179  // that we are using 'int' as the key here rather than
180  // subdomain_id_type or boundary_id_type, since at this point, it
181  // could be either.
182  typedef std::pair<unsigned, std::string> GmshPhysical;
183  std::map<int, GmshPhysical> gmsh_physicals;
184 
185  // map to hold the node numbers for translation
186  // note the the nodes can be non-consecutive
187  std::map<unsigned int, unsigned int> nodetrans;
188 
189  // Map from entity tag to physical id. The key is a pair with the first
190  // item being the dimension of the entity and the second item being
191  // the entity tag/id
192  std::map<std::pair<unsigned, int>, int> entity_to_physical_id;
193 
194  // For reading the file line by line
195  std::string s;
196 
197  while (true)
198  {
199  // Try to read something. This may set EOF!
200  std::getline(in, s);
201 
202  if (in)
203  {
204  // Process s...
205 
206  if (s.find("$MeshFormat") == static_cast<std::string::size_type>(0))
207  {
208  in >> version >> format >> size;
209  if (version < 2.0)
210  {
211  // Some notes on gmsh mesh versions:
212  //
213  // Mesh version 2.0 goes back as far as I know. It's not explicitly
214  // mentioned here: http://www.geuz.org/gmsh/doc/VERSIONS.txt
215  //
216  // As of gmsh-2.4.0:
217  // bumped mesh version format to 2.1 (small change in the $PhysicalNames
218  // section, where the group dimension is now required);
219  // [Since we don't even parse the PhysicalNames section at the time
220  // of this writing, I don't think this change affects us.]
221  //
222  // Mesh version 2.2 tested by Manav Bhatia; no other
223  // libMesh code changes were required for support
224  //
225  // Mesh version 4.0 is a near complete rewrite of the previous mesh version
226  //
227  libmesh_error_msg("Error: Unknown msh file version " << version);
228  }
229 
230  if (format)
231  libmesh_error_msg("Error: Unknown data format for mesh in Gmsh reader.");
232  }
233 
234  // Read and process the "PhysicalNames" section.
235  else if (s.find("$PhysicalNames") == static_cast<std::string::size_type>(0))
236  {
237  // The lines in the PhysicalNames section should look like the following:
238  // 2 1 "frac" lower_dimensional_block
239  // 2 3 "top"
240  // 2 4 "bottom"
241  // 3 2 "volume"
242 
243  // Read in the number of physical groups to expect in the file.
244  unsigned int num_physical_groups = 0;
245  in >> num_physical_groups;
246 
247  // Read rest of line including newline character.
248  std::getline(in, s);
249 
250  for (unsigned int i=0; i<num_physical_groups; ++i)
251  {
252  // Read an entire line of the PhysicalNames section.
253  std::getline(in, s);
254 
255  // Use an istringstream to extract the physical
256  // dimension, physical id, and physical name from
257  // this line.
258  std::istringstream s_stream(s);
259  unsigned phys_dim;
260  int phys_id;
261  std::string phys_name;
262  s_stream >> phys_dim >> phys_id >> phys_name;
263 
264  // Not sure if this is true for all Gmsh files, but
265  // my test file has quotes around the phys_name
266  // string. So let's erase any quotes now...
267  phys_name.erase(std::remove(phys_name.begin(), phys_name.end(), '"'), phys_name.end());
268 
269  // Record this ID for later assignment of subdomain/sideset names.
270  gmsh_physicals[phys_id] = std::make_pair(phys_dim, phys_name);
271 
272  // If 's' also contains the libmesh-specific string
273  // "lower_dimensional_block", add this block ID to
274  // the list of blocks which are not boundary
275  // conditions.
276  if (s.find("lower_dimensional_block") != std::string::npos)
277  {
278  lower_dimensional_blocks.insert(cast_int<subdomain_id_type>(phys_id));
279 
280  // The user has explicitly told us that this
281  // block is a subdomain, so set that association
282  // in the Mesh.
283  mesh.subdomain_name(cast_int<subdomain_id_type>(phys_id)) = phys_name;
284  }
285  }
286  }
287 
288  else if (s.find("$Entities") == static_cast<std::string::size_type>(0))
289  {
290  if (version >= 4.0)
291  {
292  std::size_t num_point_entities, num_curve_entities, num_surface_entities, num_volume_entities;
293  in >> num_point_entities >> num_curve_entities >> num_surface_entities >> num_volume_entities;
294 
295  for (std::size_t n = 0; n < num_point_entities; ++n)
296  {
297  int point_tag, physical_tag;
298  Real x, y, z;
299  std::size_t num_physical_tags;
300  in >> point_tag >> x >> y >> z >> num_physical_tags;
301  if (num_physical_tags > 1)
302  libmesh_error_msg("Sorry, you cannot currently specify multiple subdomain or " <<
303  "boundary ids for a given geometric entity");
304  else if (num_physical_tags)
305  {
306  in >> physical_tag;
307  entity_to_physical_id[std::make_pair(0, point_tag)] = physical_tag;
308  }
309  }
310  for (std::size_t n = 0; n < num_curve_entities; ++n)
311  {
312  int curve_tag, physical_tag;
313  Real minx, miny, minz, maxx, maxy, maxz;
314  std::size_t num_physical_tags;
315  in >> curve_tag >> minx >> miny >> minz >> maxx >> maxy >> maxz >> num_physical_tags;
316  if (num_physical_tags > 1)
317  libmesh_error_msg("I don't believe that we can specify multiple subdomain or " <<
318  "boundary ids for a given geometric entity");
319  else if (num_physical_tags)
320  {
321  in >> physical_tag;
322  entity_to_physical_id[std::make_pair(1, curve_tag)] = physical_tag;
323  }
324 
325  // Read to end of line; this captures bounding information that we don't care about
326  std::getline(in, s);
327  }
328  for (std::size_t n = 0; n < num_surface_entities; ++n)
329  {
330  int surface_tag, physical_tag;
331  Real minx, miny, minz, maxx, maxy, maxz;
332  std::size_t num_physical_tags;
333  in >> surface_tag >> minx >> miny >> minz >> maxx >> maxy >> maxz >> num_physical_tags;
334  if (num_physical_tags > 1)
335  libmesh_error_msg("I don't believe that we can specify multiple subdomain or " <<
336  "boundary ids for a given geometric entity");
337  else if (num_physical_tags)
338  {
339  in >> physical_tag;
340  entity_to_physical_id[std::make_pair(2, surface_tag)] = physical_tag;
341  }
342 
343  // Read to end of line; this captures bounding information that we don't care about
344  std::getline(in, s);
345  }
346  for (std::size_t n = 0; n < num_volume_entities; ++n)
347  {
348  int volume_tag, physical_tag;
349  Real minx, miny, minz, maxx, maxy, maxz;
350  std::size_t num_physical_tags;
351  in >> volume_tag >> minx >> miny >> minz >> maxx >> maxy >> maxz >> num_physical_tags;
352  if (num_physical_tags > 1)
353  libmesh_error_msg("I don't believe that we can specify multiple subdomain or " <<
354  "boundary ids for a given geometric entity");
355  else if (num_physical_tags)
356  {
357  in >> physical_tag;
358  entity_to_physical_id[std::make_pair(3, volume_tag)] = physical_tag;
359  }
360 
361  // Read to end of line; this captures bounding information that we don't care about
362  std::getline(in, s);
363  }
364  // Read the $EndEntities
365  std::getline(in, s);
366  } // end if (version >= 4.0)
367 
368  else
369  libmesh_error_msg("The Entities block was introduced in mesh version 4.0");
370  } // end if $Entities
371 
372  // read the node block
373  else if (s.find("$NOD") == static_cast<std::string::size_type>(0) ||
374  s.find("$NOE") == static_cast<std::string::size_type>(0) ||
375  s.find("$Nodes") == static_cast<std::string::size_type>(0))
376  {
377  if (version < 4.0)
378  {
379  unsigned int num_nodes = 0;
380  in >> num_nodes;
381  mesh.reserve_nodes (num_nodes);
382 
383  // read in the nodal coordinates and form points.
384  Real x, y, z;
385  unsigned int id;
386 
387  // add the nodal coordinates to the mesh
388  for (unsigned int i=0; i<num_nodes; ++i)
389  {
390  in >> id >> x >> y >> z;
391  mesh.add_point (Point(x, y, z), i);
392  nodetrans[id] = i;
393  }
394  }
395  else
396  {
397  // Read numEntityBlocks line
398  std::size_t num_entities = 0, num_nodes = 0, min_node_tag, max_node_tag;
399  in >> num_entities >> num_nodes >> min_node_tag >> max_node_tag;
400 
401  mesh.reserve_nodes(num_nodes);
402 
403  std::size_t node_counter = 0;
404 
405  // Now loop over entities
406  for (std::size_t i = 0; i < num_entities; ++i)
407  {
408  int entity_dim, entity_tag, parametric;
409  std::size_t num_nodes_in_block = 0;
410  in >> entity_dim >> entity_tag >> parametric >> num_nodes_in_block;
411  if (parametric)
412  libmesh_error_msg("We don't currently support reading parametric gmsh entities");
413 
414  // Read the node tags/ids
415  std::size_t gmsh_id;
416  for (std::size_t n = 0; n < num_nodes_in_block; ++n)
417  {
418 
419  in >> gmsh_id;
420  nodetrans[gmsh_id] = node_counter++;
421  }
422 
423  // Read the node coordinates and add the nodes to the mesh
424  Real x, y, z;
425  for (std::size_t libmesh_id = node_counter - num_nodes_in_block;
426  libmesh_id < node_counter;
427  ++libmesh_id)
428  {
429  in >> x >> y >> z;
430  mesh.add_point(Point(x, y, z), libmesh_id);
431  }
432  }
433  }
434  // read the $ENDNOD delimiter
435  std::getline(in, s);
436  }
437 
438  // Read the element block
439  else if (s.find("$ELM") == static_cast<std::string::size_type>(0) ||
440  s.find("$Elements") == static_cast<std::string::size_type>(0))
441  {
442  // Keep track of element dimensions seen
443  std::vector<unsigned> elem_dimensions_seen(3);
444 
445  if (version < 4.0)
446  {
447  // For reading the number of elements and the node ids from the stream
448  unsigned int
449  num_elem = 0,
450  node_id = 0;
451 
452  // read how many elements are there, and reserve space in the mesh
453  in >> num_elem;
454  mesh.reserve_elem (num_elem);
455 
456  // As of version 2.2, the format for each element line is:
457  // elm-number elm-type number-of-tags < tag > ... node-number-list
458  // From the Gmsh docs:
459  // * the first tag is the number of the
460  // physical entity to which the element belongs
461  // * the second is the number of the elementary geometrical
462  // entity to which the element belongs
463  // * the third is the number of mesh partitions to which the element
464  // belongs
465  // * The rest of the tags are the partition ids (negative
466  // partition ids indicate ghost cells). A zero tag is
467  // equivalent to no tag. Gmsh and most codes using the
468  // MSH 2 format require at least the first two tags
469  // (physical and elementary tags).
470 
471  // read the elements
472  for (unsigned int iel=0; iel<num_elem; ++iel)
473  {
474  unsigned int
475  id, type,
476  physical=1, elementary=1,
477  nnodes=0, ntags;
478 
479  // Note: tag has to be an int because it could be negative,
480  // see above.
481  int tag;
482 
483  if (version <= 1.0)
484  in >> id >> type >> physical >> elementary >> nnodes;
485 
486  else
487  {
488  in >> id >> type >> ntags;
489 
490  if (ntags > 2)
491  libmesh_do_once(libMesh::err << "Warning, ntags=" << ntags << ", but we currently only support reading 2 flags." << std::endl;);
492 
493  for (unsigned int j = 0; j < ntags; j++)
494  {
495  in >> tag;
496  if (j == 0)
497  physical = tag;
498  else if (j == 1)
499  elementary = tag;
500  }
501  }
502 
503  // Get a reference to the ElementDefinition, throw an error if not found.
504  const GmshIO::ElementDefinition & eletype =
505  libmesh_map_find(_element_maps.in, type);
506 
507  // If we read nnodes, make sure it matches the number in eletype.nnodes
508  if (nnodes != 0 && nnodes != eletype.nnodes)
509  libmesh_error_msg("nnodes = " << nnodes << " and eletype.nnodes = " << eletype.nnodes << " do not match.");
510 
511  // Assign the value from the eletype object.
512  nnodes = eletype.nnodes;
513 
514  // Don't add 0-dimensional "point" elements to the
515  // Mesh. They should *always* be treated as boundary
516  // "nodeset" data.
517  if (eletype.dim > 0)
518  {
519  // Record this element dimension as being "seen".
520  // We will treat all elements with dimension <
521  // max(dimension) as specifying boundary conditions,
522  // but we won't know what max_elem_dimension_seen is
523  // until we read the entire file.
524  elem_dimensions_seen[eletype.dim-1] = 1;
525 
526  // Add the element to the mesh
527  {
528  Elem * elem = Elem::build(eletype.type).release();
529  elem->set_id(iel);
530  elem = mesh.add_elem(elem);
531 
532  // Make sure that the libmesh element we added has nnodes nodes.
533  if (elem->n_nodes() != nnodes)
534  libmesh_error_msg("Number of nodes for element " \
535  << id \
536  << " of type " << eletype.type \
537  << " (Gmsh type " << type \
538  << ") does not match Libmesh definition. " \
539  << "I expected " << elem->n_nodes() \
540  << " nodes, but got " << nnodes);
541 
542  // Add node pointers to the elements.
543  // If there is a node translation table, use it.
544  if (eletype.nodes.size() > 0)
545  for (unsigned int i=0; i<nnodes; i++)
546  {
547  in >> node_id;
548  elem->set_node(eletype.nodes[i]) = mesh.node_ptr(nodetrans[node_id]);
549  }
550  else
551  {
552  for (unsigned int i=0; i<nnodes; i++)
553  {
554  in >> node_id;
555  elem->set_node(i) = mesh.node_ptr(nodetrans[node_id]);
556  }
557  }
558 
559  // Finally, set the subdomain ID to physical. If this is a lower-dimension element, this ID will
560  // eventually go into the Mesh's BoundaryInfo object.
561  elem->subdomain_id() = static_cast<subdomain_id_type>(physical);
562  }
563  }
564 
565  // Handle 0-dimensional elements (points) by adding
566  // them to the BoundaryInfo object with
567  // boundary_id == physical.
568  else
569  {
570  // This seems like it should always be the same
571  // number as the 'id' we already read in on this
572  // line. At least it was in the example gmsh
573  // file I had...
574  in >> node_id;
576  (nodetrans[node_id],
577  static_cast<boundary_id_type>(physical));
578  }
579  } // element loop
580  } // end if (version < 4.0)
581 
582  else
583  {
584  std::size_t num_entity_blocks = 0, num_elem = 0, min_element_tag, max_element_tag;
585 
586  // Read entity information
587  in >> num_entity_blocks >> num_elem >> min_element_tag >> max_element_tag;
588 
589  mesh.reserve_elem(num_elem);
590 
591  std::size_t iel = 0;
592 
593  // Loop over entity blocks
594  for (std::size_t i = 0; i < num_entity_blocks; ++i)
595  {
596  int entity_dim, entity_tag;
597  unsigned int element_type;
598  std::size_t num_elems_in_block = 0;
599  in >> entity_dim >> entity_tag >> element_type >> num_elems_in_block;
600 
601  // Get a reference to the ElementDefinition
602  const GmshIO::ElementDefinition & eletype =
603  libmesh_map_find(_element_maps.in, element_type);
604 
605  // Don't add 0-dimensional "point" elements to the
606  // Mesh. They should *always* be treated as boundary
607  // "nodeset" data.
608  if (eletype.dim > 0)
609  {
610  // Record this element dimension as being "seen".
611  // We will treat all elements with dimension <
612  // max(dimension) as specifying boundary conditions,
613  // but we won't know what max_elem_dimension_seen is
614  // until we read the entire file.
615  elem_dimensions_seen[eletype.dim-1] = 1;
616 
617  // Loop over elements with dim > 0
618  for (std::size_t n = 0; n < num_elems_in_block; ++n)
619  {
620  Elem * elem = Elem::build(eletype.type).release();
621  elem->set_id(iel++);
622  elem = mesh.add_elem(elem);
623 
624  std::size_t gmsh_element_id;
625  in >> gmsh_element_id;
626 
627  // Get the remainer of the line that includes the nodes ids
628  std::getline(in, s);
629  std::istringstream is(s);
630  std::size_t local_node_counter = 0, gmsh_node_id;
631  while (is >> gmsh_node_id)
632  {
633  // Add node pointers to the elements.
634  // If there is a node translation table, use it.
635  if (eletype.nodes.size() > 0)
636  elem->set_node(eletype.nodes[local_node_counter++]) =
637  mesh.node_ptr(nodetrans[gmsh_node_id]);
638  else
639  elem->set_node(local_node_counter++) = mesh.node_ptr(nodetrans[gmsh_node_id]);
640  }
641 
642  // Make sure that the libmesh element we added has nnodes nodes.
643  if (elem->n_nodes() != local_node_counter)
644  libmesh_error_msg("Number of nodes for element " \
645  << gmsh_element_id \
646  << " of type " << eletype.type \
647  << " (Gmsh type " << element_type \
648  << ") does not match Libmesh definition. " \
649  << "I expected " << elem->n_nodes() \
650  << " nodes, but got " << local_node_counter);
651 
652 
653  // Finally, set the subdomain ID to physical. If this is a lower-dimension element, this ID will
654  // eventually go into the Mesh's BoundaryInfo object.
655  elem->subdomain_id() = static_cast<subdomain_id_type>(
656  entity_to_physical_id[std::make_pair(entity_dim, entity_tag)]);
657 
658  } // end for (loop over elements in entity block)
659  } // end if (eletype.dim > 0)
660 
661  else
662  {
663  for (std::size_t n = 0; n < num_elems_in_block; ++n)
664  {
665  std::size_t gmsh_element_id, gmsh_node_id;
666  in >> gmsh_element_id;
667  in >> gmsh_node_id;
669  nodetrans[gmsh_node_id],
670  static_cast<boundary_id_type>(entity_to_physical_id[
671  std::make_pair(entity_dim, entity_tag)]));
672  } // end for (loop over elements in entity block)
673  } // end if (eletype.dim == 0)
674  } // end for (loop over entity blocks)
675  } // end if (version >= 4.0)
676 
677  // read the $ENDELM delimiter
678  std::getline(in, s);
679 
680  // Record the max and min element dimension seen while reading the file.
681  unsigned char
682  max_elem_dimension_seen=1,
683  min_elem_dimension_seen=3;
684 
685  for (auto i : index_range(elem_dimensions_seen))
686  if (elem_dimensions_seen[i])
687  {
688  // Debugging
689  // libMesh::out << "Seen elements of dimension " << i+1 << std::endl;
690  max_elem_dimension_seen =
691  std::max(max_elem_dimension_seen, cast_int<unsigned char>(i+1));
692  min_elem_dimension_seen =
693  std::min(min_elem_dimension_seen, cast_int<unsigned char>(i+1));
694  }
695 
696  // Debugging:
697  // libMesh::out << "max_elem_dimension_seen=" << max_elem_dimension_seen << std::endl;
698  // libMesh::out << "min_elem_dimension_seen=" << min_elem_dimension_seen << std::endl;
699 
700 
701  // How many different element dimensions did we see while reading from file?
702  unsigned n_dims_seen = std::accumulate(elem_dimensions_seen.begin(),
703  elem_dimensions_seen.end(),
704  static_cast<unsigned>(0),
705  std::plus<unsigned>());
706 
707 
708  // Set mesh_dimension based on the largest element dimension seen.
709  mesh.set_mesh_dimension(max_elem_dimension_seen);
710 
711  // Now that we know the maximum element dimension seen,
712  // we know whether the physical names are subdomain
713  // names or sideset names.
714  for (const auto & pr : gmsh_physicals)
715  {
716  // Extract data
717  int phys_id = pr.first;
718  unsigned phys_dim = pr.second.first;
719  const std::string & phys_name = pr.second.second;
720 
721  // If the physical's dimension matches the largest
722  // dimension we've seen, it's a subdomain name.
723  if (phys_dim == max_elem_dimension_seen)
724  mesh.subdomain_name(cast_int<subdomain_id_type>(phys_id)) = phys_name;
725 
726  // Otherwise, if it's not a lower-dimensional
727  // block, it's a sideset name.
728  else if (phys_dim < max_elem_dimension_seen &&
729  !lower_dimensional_blocks.count(cast_int<boundary_id_type>(phys_id)))
730  mesh.get_boundary_info().sideset_name(cast_int<boundary_id_type>(phys_id)) = phys_name;
731  }
732 
733  if (n_dims_seen > 1)
734  {
735  // Store lower-dimensional elements in a map sorted
736  // by Elem::key(). We use a multimap for two reasons:
737  // 1.) The hash function is not guaranteed to be
738  // unique, so different lower-dimensional elements
739  // could theoretically hash to the same value,
740  // although this is pretty unlikely.
741  // 2.) The Gmsh file may contain multiple
742  // lower-dimensional elements for a single side in
743  // order to implement multiple boundary ids for a
744  // single side. These lower-dimensional elements
745  // will all hash to the same value, and we need to
746  // be able to store all of them.
747  typedef std::unordered_multimap<dof_id_type, Elem *> provide_container_t;
748  provide_container_t provide_bcs;
749 
750  // 1st loop over active elements - get info about lower-dimensional elements.
751  for (auto & elem : mesh.active_element_ptr_range())
752  if (elem->dim() < max_elem_dimension_seen &&
753  !lower_dimensional_blocks.count(elem->subdomain_id()))
754  {
755  // To be consistent with the previous
756  // GmshIO behavior, add all the
757  // lower-dimensional elements' nodes to
758  // the Mesh's BoundaryInfo object with the
759  // lower-dimensional element's subdomain
760  // ID.
761  for (auto n : elem->node_index_range())
762  mesh.get_boundary_info().add_node(elem->node_id(n),
763  elem->subdomain_id());
764 
765  // Store this elem in a quickly-searchable
766  // container to use it to assign boundary
767  // conditions later.
768  provide_bcs.insert(std::make_pair(elem->key(), elem));
769  }
770 
771  // 2nd loop over active elements - use lower dimensional element data to set BCs for higher dimensional elements
772  for (auto & elem : mesh.active_element_ptr_range())
773  if (elem->dim() == max_elem_dimension_seen)
774  {
775  // This is a max-dimension element that
776  // may require BCs. For each of its
777  // sides, including internal sides, we'll
778  // see if one more more lower-dimensional elements
779  // provides boundary information for it.
780  // Note that we have not yet called
781  // find_neighbors(), so we can't use
782  // elem->neighbor(sn) in this algorithm...
783  for (auto sn : elem->side_index_range())
784  for (const auto & pr : as_range(provide_bcs.equal_range(elem->key(sn))))
785  {
786  // For each side side in the provide_bcs multimap...
787  // Construct the side for hash verification.
788  std::unique_ptr<Elem> side (elem->build_side_ptr(sn));
789 
790  // Construct the lower-dimensional element to compare to the side.
791  Elem * lower_dim_elem = pr.second;
792 
793  // This was a hash, so it might not be perfect. Let's verify...
794  if (*lower_dim_elem == *side)
795  {
796  // Add the lower-dimensional
797  // element's subdomain_id as a
798  // boundary_id for the
799  // higher-dimensional element.
800  boundary_id_type bid = cast_int<boundary_id_type>(lower_dim_elem->subdomain_id());
801  mesh.get_boundary_info().add_side(elem, sn, bid);
802  }
803  }
804  }
805 
806  // 3rd loop over active elements - Remove the lower-dimensional elements
807  for (auto & elem : mesh.active_element_ptr_range())
808  if (elem->dim() < max_elem_dimension_seen &&
809  !lower_dimensional_blocks.count(elem->subdomain_id()))
810  mesh.delete_elem(elem);
811  } // end if (n_dims_seen > 1)
812  } // if $ELM
813 
814  continue;
815  } // if (in)
816 
817 
818  // If !in, check to see if EOF was set. If so, break out
819  // of while loop.
820  if (in.eof())
821  break;
822 
823  // If !in and !in.eof(), stream is in a bad state!
824  libmesh_error_msg("Stream is bad! Perhaps the file does not exist?");
825 
826  } // while true
827 }
828 
829 
830 
831 void GmshIO::write (const std::string & name)
832 {
833  if (MeshOutput<MeshBase>::mesh().processor_id() == 0)
834  {
835  // Open the output file stream
836  std::ofstream out_stream (name.c_str());
837 
838  // Make sure it opened correctly
839  if (!out_stream.good())
840  libmesh_file_error(name.c_str());
841 
842  this->write_mesh (out_stream);
843  }
844 }
845 
846 
847 
848 void GmshIO::write_nodal_data (const std::string & fname,
849  const std::vector<Number> & soln,
850  const std::vector<std::string> & names)
851 {
852  LOG_SCOPE("write_nodal_data()", "GmshIO");
853 
854  if (MeshOutput<MeshBase>::mesh().processor_id() == 0)
855  this->write_post (fname, &soln, &names);
856 }
857 
858 
859 
860 void GmshIO::write_mesh (std::ostream & out_stream)
861 {
862  // Be sure that the stream is valid.
863  libmesh_assert (out_stream.good());
864 
865  // Get a const reference to the mesh
867 
868  // If requested, write out lower-dimensional elements for
869  // element-side-based boundary conditions.
870  std::size_t n_boundary_faces = 0;
872  n_boundary_faces = mesh.get_boundary_info().n_boundary_conds();
873 
874  // Note: we are using version 2.0 of the gmsh output format.
875 
876  // Write the file header.
877  out_stream << "$MeshFormat\n";
878  out_stream << "2.0 0 " << sizeof(Real) << '\n';
879  out_stream << "$EndMeshFormat\n";
880 
881  // write the nodes in (n x y z) format
882  out_stream << "$Nodes\n";
883  out_stream << mesh.n_nodes() << '\n';
884 
885  for (auto v : IntRange<unsigned int>(0, mesh.n_nodes()))
886  out_stream << mesh.node_ref(v).id()+1 << " "
887  << mesh.node_ref(v)(0) << " "
888  << mesh.node_ref(v)(1) << " "
889  << mesh.node_ref(v)(2) << '\n';
890  out_stream << "$EndNodes\n";
891 
892  {
893  // write the connectivity
894  out_stream << "$Elements\n";
895  out_stream << mesh.n_active_elem() + n_boundary_faces << '\n';
896 
897  // loop over the elements
898  for (const auto & elem : mesh.active_element_ptr_range())
899  {
900  // Get a reference to the ElementDefinition object
901  const ElementDefinition & eletype =
902  libmesh_map_find(_element_maps.out, elem->type());
903 
904  // The element mapper better not require any more nodes
905  // than are present in the current element!
906  libmesh_assert_less_equal (eletype.nodes.size(), elem->n_nodes());
907 
908  // elements ids are 1 based in Gmsh
909  out_stream << elem->id()+1 << " ";
910  // element type
911  out_stream << eletype.gmsh_type;
912 
913  // write the number of tags (3) and their values:
914  // 1 (physical entity)
915  // 2 (geometric entity)
916  // 3 (partition entity)
917  out_stream << " 3 "
918  << static_cast<unsigned int>(elem->subdomain_id())
919  << " 0 "
920  << elem->processor_id()+1
921  << " ";
922 
923  // if there is a node translation table, use it
924  if (eletype.nodes.size() > 0)
925  for (auto i : elem->node_index_range())
926  out_stream << elem->node_id(eletype.nodes[i])+1 << " "; // gmsh is 1-based
927  // otherwise keep the same node order
928  else
929  for (auto i : elem->node_index_range())
930  out_stream << elem->node_id(i)+1 << " "; // gmsh is 1-based
931  out_stream << "\n";
932  } // element loop
933  }
934 
935  {
936  // A counter for writing surface elements to the Gmsh file
937  // sequentially. We start numbering them with a number strictly
938  // larger than the largest element ID in the mesh. Note: the
939  // MeshBase docs say "greater than or equal to" the maximum
940  // element id in the mesh, so technically we might need a +1 here,
941  // but all of the implementations return an ID strictly greater
942  // than the largest element ID in the Mesh.
943  unsigned int e_id = mesh.max_elem_id();
944 
945  // loop over the elements, writing out boundary faces
946  if (n_boundary_faces)
947  {
948  // Build a list of (elem, side, bc) tuples.
949  auto bc_triples = mesh.get_boundary_info().build_side_list();
950 
951  // Loop over these lists, writing data to the file.
952  for (const auto & t : bc_triples)
953  {
954  const Elem & elem = mesh.elem_ref(std::get<0>(t));
955 
956  std::unique_ptr<const Elem> side = elem.build_side_ptr(std::get<1>(t));
957 
958  // consult the export element table
959  const GmshIO::ElementDefinition & eletype =
960  libmesh_map_find(_element_maps.out, side->type());
961 
962  // The element mapper better not require any more nodes
963  // than are present in the current element!
964  libmesh_assert_less_equal (eletype.nodes.size(), side->n_nodes());
965 
966  // elements ids are 1-based in Gmsh
967  out_stream << e_id+1 << " ";
968 
969  // element type
970  out_stream << eletype.gmsh_type;
971 
972  // write the number of tags:
973  // 1 (physical entity)
974  // 2 (geometric entity)
975  // 3 (partition entity)
976  out_stream << " 3 "
977  << std::get<2>(t)
978  << " 0 "
979  << elem.processor_id()+1
980  << " ";
981 
982  // if there is a node translation table, use it
983  if (eletype.nodes.size() > 0)
984  for (auto i : side->node_index_range())
985  out_stream << side->node_id(eletype.nodes[i])+1 << " "; // gmsh is 1-based
986 
987  // otherwise keep the same node order
988  else
989  for (auto i : side->node_index_range())
990  out_stream << side->node_id(i)+1 << " "; // gmsh is 1-based
991 
992  // Go to the next line
993  out_stream << "\n";
994 
995  // increment this index too...
996  ++e_id;
997  }
998  }
999 
1000  out_stream << "$EndElements\n";
1001  }
1002 }
1003 
1004 
1005 
1006 void GmshIO::write_post (const std::string & fname,
1007  const std::vector<Number> * v,
1008  const std::vector<std::string> * solution_names)
1009 {
1010 
1011  // Should only do this on processor 0!
1012  libmesh_assert_equal_to (MeshOutput<MeshBase>::mesh().processor_id(), 0);
1013 
1014  // Create an output stream
1015  std::ofstream out_stream(fname.c_str());
1016 
1017  // Make sure it opened correctly
1018  if (!out_stream.good())
1019  libmesh_file_error(fname.c_str());
1020 
1021  // create a character buffer
1022  char buf[80];
1023 
1024  // Get a constant reference to the mesh.
1026 
1027  // write the data
1028  if ((solution_names != nullptr) && (v != nullptr))
1029  {
1030  const unsigned int n_vars =
1031  cast_int<unsigned int>(solution_names->size());
1032 
1033  if (!(v->size() == mesh.n_nodes()*n_vars))
1034  libMesh::err << "ERROR: v->size()=" << v->size()
1035  << ", mesh.n_nodes()=" << mesh.n_nodes()
1036  << ", n_vars=" << n_vars
1037  << ", mesh.n_nodes()*n_vars=" << mesh.n_nodes()*n_vars
1038  << "\n";
1039 
1040  libmesh_assert_equal_to (v->size(), mesh.n_nodes()*n_vars);
1041 
1042  // write the header
1043  out_stream << "$PostFormat\n";
1044  if (this->binary())
1045  out_stream << "1.2 1 " << sizeof(double) << "\n";
1046  else
1047  out_stream << "1.2 0 " << sizeof(double) << "\n";
1048  out_stream << "$EndPostFormat\n";
1049 
1050  // Loop over the elements to see how much of each type there are
1051  unsigned int n_points=0, n_lines=0, n_triangles=0, n_quadrangles=0,
1052  n_tetrahedra=0, n_hexahedra=0, n_prisms=0, n_pyramids=0;
1053  unsigned int n_scalar=0, n_vector=0, n_tensor=0;
1054  unsigned int nb_text2d=0, nb_text2d_chars=0, nb_text3d=0, nb_text3d_chars=0;
1055 
1056  {
1057  for (const auto & elem : mesh.active_element_ptr_range())
1058  {
1059  const ElemType elemtype = elem->type();
1060 
1061  switch (elemtype)
1062  {
1063  case EDGE2:
1064  case EDGE3:
1065  case EDGE4:
1066  {
1067  n_lines += 1;
1068  break;
1069  }
1070  case TRI3:
1071  case TRI6:
1072  {
1073  n_triangles += 1;
1074  break;
1075  }
1076  case QUAD4:
1077  case QUAD8:
1078  case QUAD9:
1079  {
1080  n_quadrangles += 1;
1081  break;
1082  }
1083  case TET4:
1084  case TET10:
1085  {
1086  n_tetrahedra += 1;
1087  break;
1088  }
1089  case HEX8:
1090  case HEX20:
1091  case HEX27:
1092  {
1093  n_hexahedra += 1;
1094  break;
1095  }
1096  case PRISM6:
1097  case PRISM15:
1098  case PRISM18:
1099  {
1100  n_prisms += 1;
1101  break;
1102  }
1103  case PYRAMID5:
1104  {
1105  n_pyramids += 1;
1106  break;
1107  }
1108  default:
1109  libmesh_error_msg("ERROR: Nonexistent element type " << elem->type());
1110  }
1111  }
1112  }
1113 
1114  // create a view for each variable
1115  for (unsigned int ivar=0; ivar < n_vars; ivar++)
1116  {
1117  std::string varname = (*solution_names)[ivar];
1118 
1119  // at the moment, we just write out scalar quantities
1120  // later this should be made configurable through
1121  // options to the writer class
1122  n_scalar = 1;
1123 
1124  // write the variable as a view, and the number of time steps
1125  out_stream << "$View\n" << varname << " " << 1 << "\n";
1126 
1127  // write how many of each geometry type are written
1128  out_stream << n_points * n_scalar << " "
1129  << n_points * n_vector << " "
1130  << n_points * n_tensor << " "
1131  << n_lines * n_scalar << " "
1132  << n_lines * n_vector << " "
1133  << n_lines * n_tensor << " "
1134  << n_triangles * n_scalar << " "
1135  << n_triangles * n_vector << " "
1136  << n_triangles * n_tensor << " "
1137  << n_quadrangles * n_scalar << " "
1138  << n_quadrangles * n_vector << " "
1139  << n_quadrangles * n_tensor << " "
1140  << n_tetrahedra * n_scalar << " "
1141  << n_tetrahedra * n_vector << " "
1142  << n_tetrahedra * n_tensor << " "
1143  << n_hexahedra * n_scalar << " "
1144  << n_hexahedra * n_vector << " "
1145  << n_hexahedra * n_tensor << " "
1146  << n_prisms * n_scalar << " "
1147  << n_prisms * n_vector << " "
1148  << n_prisms * n_tensor << " "
1149  << n_pyramids * n_scalar << " "
1150  << n_pyramids * n_vector << " "
1151  << n_pyramids * n_tensor << " "
1152  << nb_text2d << " "
1153  << nb_text2d_chars << " "
1154  << nb_text3d << " "
1155  << nb_text3d_chars << "\n";
1156 
1157  // if binary, write a marker to identify the endianness of the file
1158  if (this->binary())
1159  {
1160  const int one = 1;
1161  std::memcpy(buf, &one, sizeof(int));
1162  out_stream.write(buf, sizeof(int));
1163  }
1164 
1165  // the time steps (there is just 1 at the moment)
1166  if (this->binary())
1167  {
1168  double one = 1;
1169  std::memcpy(buf, &one, sizeof(double));
1170  out_stream.write(buf, sizeof(double));
1171  }
1172  else
1173  out_stream << "1\n";
1174 
1175  // Loop over the elements and write out the data
1176  for (const auto & elem : mesh.active_element_ptr_range())
1177  {
1178  const unsigned int nv = elem->n_vertices();
1179  // this is quite crappy, but I did not invent that file format!
1180  for (unsigned int d=0; d<3; d++) // loop over the dimensions
1181  {
1182  for (unsigned int n=0; n < nv; n++) // loop over vertices
1183  {
1184  const Point & vertex = elem->point(n);
1185  if (this->binary())
1186  {
1187 #if defined(LIBMESH_DEFAULT_TRIPLE_PRECISION) || defined(LIBMESH_DEFAULT_QUADRUPLE_PRECISION)
1188  libmesh_warning("Gmsh binary writes use only double precision!");
1189 #endif
1190  double tmp = double(vertex(d));
1191  std::memcpy(buf, &tmp, sizeof(double));
1192  out_stream.write(reinterpret_cast<char *>(buf), sizeof(double));
1193  }
1194  else
1195  out_stream << vertex(d) << " ";
1196  }
1197  if (!this->binary())
1198  out_stream << "\n";
1199  }
1200 
1201  // now finally write out the data
1202  for (unsigned int i=0; i < nv; i++) // loop over vertices
1203  if (this->binary())
1204  {
1205 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1206  libMesh::out << "WARNING: Gmsh::write_post does not fully support "
1207  << "complex numbers. Will only write the real part of "
1208  << "variable " << varname << std::endl;
1209 #endif
1210  double tmp = double(libmesh_real((*v)[elem->node_id(i)*n_vars + ivar]));
1211  std::memcpy(buf, &tmp, sizeof(double));
1212  out_stream.write(reinterpret_cast<char *>(buf), sizeof(double));
1213  }
1214  else
1215  {
1216 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1217  libMesh::out << "WARNING: Gmsh::write_post does not fully support "
1218  << "complex numbers. Will only write the real part of "
1219  << "variable " << varname << std::endl;
1220 #endif
1221  out_stream << libmesh_real((*v)[elem->node_id(i)*n_vars + ivar]) << "\n";
1222  }
1223  }
1224  if (this->binary())
1225  out_stream << "\n";
1226  out_stream << "$EndView\n";
1227 
1228  } // end variable loop (writing the views)
1229  }
1230 }
1231 
1232 } // namespace libMesh
libMesh::GmshIO::ElementMaps::in
std::map< unsigned int, ElementDefinition > in
Definition: gmsh_io.h:193
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::GmshIO::read
virtual void read(const std::string &name) override
Reads in a mesh in the Gmsh *.msh format from the ASCII file given by name.
Definition: gmsh_io.C:147
libMesh::MeshBase::reserve_nodes
virtual void reserve_nodes(const dof_id_type nn)=0
Reserves space for a known number of nodes.
libMesh::GmshIO::ElementDefinition::dim
unsigned int dim
Definition: gmsh_io.h:174
libMesh::Elem::build_side_ptr
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy=true)=0
libMesh::GmshIO::_binary
bool _binary
Flag to write binary data.
Definition: gmsh_io.h:149
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
n_vars
unsigned int n_vars
Definition: adaptivity_ex3.C:116
libMesh::GmshIO::ElementMaps
struct which holds a map from Gmsh to libMesh element numberings and vice-versa.
Definition: gmsh_io.h:183
libMesh::Elem::n_nodes
virtual unsigned int n_nodes() const =0
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::BoundaryInfo::add_node
void add_node(const Node *node, const boundary_id_type id)
Add Node node with boundary id id to the boundary information data structures.
Definition: boundary_info.C:636
libMesh::EDGE4
Definition: enum_elem_type.h:37
libMesh::BoundaryInfo::n_boundary_conds
std::size_t n_boundary_conds() const
Definition: boundary_info.C:1615
libMesh::libmesh_real
T libmesh_real(T a)
Definition: libmesh_common.h:166
libMesh::MeshBase::delete_elem
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
libMesh::GmshIO::ElementDefinition::nnodes
unsigned int nnodes
Definition: gmsh_io.h:175
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::MeshBase::elem_ref
virtual const Elem & elem_ref(const dof_id_type i) const
Definition: mesh_base.h:521
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh::MeshBase::max_elem_id
virtual dof_id_type max_elem_id() const =0
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::GmshIO::ElementDefinition::nodes
std::vector< unsigned int > nodes
Definition: gmsh_io.h:176
libMesh::TET10
Definition: enum_elem_type.h:46
libMesh::GmshIO::_element_maps
static ElementMaps _element_maps
A static ElementMaps object that is built statically and used by all instances of this class.
Definition: gmsh_io.h:200
libMesh::GmshIO::write_lower_dimensional_elements
bool & write_lower_dimensional_elements()
Access to the flag which controls whether boundary elements are written to the Mesh file.
Definition: gmsh_io.C:140
libMesh::BoundaryInfo::sideset_name
std::string & sideset_name(boundary_id_type id)
Definition: boundary_info.C:2341
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::GmshIO::_write_lower_dimensional_elements
bool _write_lower_dimensional_elements
If true, lower-dimensional elements based on the boundary conditions get written to the output file.
Definition: gmsh_io.h:155
libMesh::GmshIO::ElementDefinition::type
ElemType type
Definition: gmsh_io.h:172
libMesh::GmshIO::write_nodal_data
virtual void write_nodal_data(const std::string &, const std::vector< Number > &, const std::vector< std::string > &) override
This method implements writing a mesh with nodal data to a specified file where the nodal data and va...
Definition: gmsh_io.C:848
libMesh::boundary_id_type
int8_t boundary_id_type
Definition: id_types.h:51
libMesh::DofObject::processor_id
processor_id_type processor_id() const
Definition: dof_object.h:829
libMesh::GmshIO::write_post
void write_post(const std::string &, const std::vector< Number > *=nullptr, const std::vector< std::string > *=nullptr)
This method implements writing a mesh with nodal data to a specified file where the nodal data and va...
Definition: gmsh_io.C:1006
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::BoundaryInfo::build_side_list
void build_side_list(std::vector< dof_id_type > &element_id_list, std::vector< unsigned short int > &side_list, std::vector< boundary_id_type > &bc_id_list) const
Creates a list of element numbers, sides, and ids for those sides.
Definition: boundary_info.C:1976
libMesh::PRISM15
Definition: enum_elem_type.h:51
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::HEX27
Definition: enum_elem_type.h:49
libMesh::GmshIO::GmshIO
GmshIO(MeshBase &mesh)
Constructor.
Definition: gmsh_io.C:123
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::GmshIO::ElementDefinition::gmsh_type
unsigned int gmsh_type
Definition: gmsh_io.h:173
libMesh::GmshIO::binary
bool & binary()
Flag indicating whether or not to write a binary file.
Definition: gmsh_io.C:133
libMesh::GmshIO::ElementMaps::add_def
void add_def(const ElementDefinition &eledef)
Definition: gmsh_io.h:186
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::GmshIO::read_mesh
void read_mesh(std::istream &in)
Implementation of the read() function.
Definition: gmsh_io.C:155
libMesh::is
PetscErrorCode PetscInt const PetscInt IS * is
Definition: petsc_dm_wrapper.C:60
libMesh::Elem::set_node
virtual Node *& set_node(const unsigned int i)
Definition: elem.h:2059
libMesh::GmshIO::write_mesh
void write_mesh(std::ostream &out)
This method implements writing a mesh to a specified file.
Definition: gmsh_io.C:860
libMesh::MeshBase::node_ref
virtual const Node & node_ref(const dof_id_type i) const
Definition: mesh_base.h:451
libMesh::MeshBase::n_nodes
virtual dof_id_type n_nodes() const =0
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::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::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::EDGE3
Definition: enum_elem_type.h:36
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::NODEELEM
Definition: enum_elem_type.h:66
libMesh::MeshInput< MeshBase >::mesh
MeshBase & 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::GmshIO::write
virtual void write(const std::string &name) override
This method implements writing a mesh to a specified file in the Gmsh *.msh format.
Definition: gmsh_io.C:831
libMesh::GmshIO::ElementMaps::out
std::map< ElemType, ElementDefinition > out
Definition: gmsh_io.h:192
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::MeshBase::clear
virtual void clear()
Deletes all the element and node data that is currently stored.
Definition: mesh_base.C:429
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::PRISM18
Definition: enum_elem_type.h:52
libMesh::GmshIO::ElementDefinition
Defines mapping from libMesh element types to Gmsh element types or vice-versa.
Definition: gmsh_io.h:160
libMesh::Elem::build
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
Definition: elem.C:246
libMesh::out
OStreamProxy out
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::MeshBase::reserve_elem
virtual void reserve_elem(const dof_id_type ne)=0
Reserves space for a known number of elements.
libMesh::GmshIO::build_element_maps
static ElementMaps build_element_maps()
A static function used to construct the _element_maps struct, statically.
Definition: gmsh_io.C:45
libMesh::MeshBase::n_active_elem
virtual dof_id_type n_active_elem() const =0
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::MeshInput
This class defines an abstract interface for Mesh input.
Definition: mesh_base.h:60
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::QUAD8
Definition: enum_elem_type.h:42
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33