23 #include <unordered_map>
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"
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);
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);
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);
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);
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);
117 _write_lower_dimensional_elements(true)
127 _write_lower_dimensional_elements(true)
149 std::ifstream in (
name.c_str());
169 int format=0, size=0;
174 std::set<subdomain_id_type> lower_dimensional_blocks;
182 typedef std::pair<unsigned, std::string> GmshPhysical;
183 std::map<int, GmshPhysical> gmsh_physicals;
187 std::map<unsigned int, unsigned int> nodetrans;
192 std::map<std::pair<unsigned, int>,
int> entity_to_physical_id;
206 if (s.find(
"$MeshFormat") == static_cast<std::string::size_type>(0))
208 in >> version >> format >> size;
227 libmesh_error_msg(
"Error: Unknown msh file version " << version);
231 libmesh_error_msg(
"Error: Unknown data format for mesh in Gmsh reader.");
235 else if (s.find(
"$PhysicalNames") == static_cast<std::string::size_type>(0))
244 unsigned int num_physical_groups = 0;
245 in >> num_physical_groups;
250 for (
unsigned int i=0; i<num_physical_groups; ++i)
258 std::istringstream s_stream(s);
261 std::string phys_name;
262 s_stream >> phys_dim >> phys_id >> phys_name;
267 phys_name.erase(std::remove(phys_name.begin(), phys_name.end(),
'"'), phys_name.end());
270 gmsh_physicals[phys_id] = std::make_pair(phys_dim, phys_name);
276 if (s.find(
"lower_dimensional_block") != std::string::npos)
278 lower_dimensional_blocks.insert(cast_int<subdomain_id_type>(phys_id));
288 else if (s.find(
"$Entities") == static_cast<std::string::size_type>(0))
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;
295 for (std::size_t n = 0; n < num_point_entities; ++n)
297 int point_tag, physical_tag;
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)
307 entity_to_physical_id[std::make_pair(0, point_tag)] = physical_tag;
310 for (std::size_t n = 0; n < num_curve_entities; ++n)
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)
322 entity_to_physical_id[std::make_pair(1, curve_tag)] = physical_tag;
328 for (std::size_t n = 0; n < num_surface_entities; ++n)
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)
340 entity_to_physical_id[std::make_pair(2, surface_tag)] = physical_tag;
346 for (std::size_t n = 0; n < num_volume_entities; ++n)
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)
358 entity_to_physical_id[std::make_pair(3, volume_tag)] = physical_tag;
369 libmesh_error_msg(
"The Entities block was introduced in mesh version 4.0");
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))
379 unsigned int num_nodes = 0;
388 for (
unsigned int i=0; i<num_nodes; ++i)
390 in >>
id >> x >> y >> z;
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;
403 std::size_t node_counter = 0;
406 for (std::size_t i = 0; i < num_entities; ++i)
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;
412 libmesh_error_msg(
"We don't currently support reading parametric gmsh entities");
416 for (std::size_t n = 0; n < num_nodes_in_block; ++n)
420 nodetrans[gmsh_id] = node_counter++;
425 for (std::size_t libmesh_id = node_counter - num_nodes_in_block;
426 libmesh_id < node_counter;
439 else if (s.find(
"$ELM") == static_cast<std::string::size_type>(0) ||
440 s.find(
"$Elements") == static_cast<std::string::size_type>(0))
443 std::vector<unsigned> elem_dimensions_seen(3);
472 for (
unsigned int iel=0; iel<num_elem; ++iel)
476 physical=1, elementary=1,
484 in >>
id >> type >> physical >> elementary >> nnodes;
488 in >>
id >> type >> ntags;
491 libmesh_do_once(
libMesh::err <<
"Warning, ntags=" << ntags <<
", but we currently only support reading 2 flags." << std::endl;);
493 for (
unsigned int j = 0; j < ntags; j++)
508 if (nnodes != 0 && nnodes != eletype.
nnodes)
509 libmesh_error_msg(
"nnodes = " << nnodes <<
" and eletype.nnodes = " << eletype.
nnodes <<
" do not match.");
524 elem_dimensions_seen[eletype.
dim-1] = 1;
534 libmesh_error_msg(
"Number of nodes for element " \
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);
544 if (eletype.
nodes.size() > 0)
545 for (
unsigned int i=0; i<nnodes; i++)
552 for (
unsigned int i=0; i<nnodes; i++)
561 elem->
subdomain_id() = static_cast<subdomain_id_type>(physical);
577 static_cast<boundary_id_type>(physical));
584 std::size_t num_entity_blocks = 0, num_elem = 0, min_element_tag, max_element_tag;
587 in >> num_entity_blocks >> num_elem >> min_element_tag >> max_element_tag;
594 for (std::size_t i = 0; i < num_entity_blocks; ++i)
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;
615 elem_dimensions_seen[eletype.
dim-1] = 1;
618 for (std::size_t n = 0; n < num_elems_in_block; ++n)
624 std::size_t gmsh_element_id;
625 in >> gmsh_element_id;
629 std::istringstream
is(s);
630 std::size_t local_node_counter = 0, gmsh_node_id;
631 while (
is >> gmsh_node_id)
635 if (eletype.
nodes.size() > 0)
643 if (elem->
n_nodes() != local_node_counter)
644 libmesh_error_msg(
"Number of nodes for element " \
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);
656 entity_to_physical_id[std::make_pair(entity_dim, entity_tag)]);
663 for (std::size_t n = 0; n < num_elems_in_block; ++n)
665 std::size_t gmsh_element_id, gmsh_node_id;
666 in >> gmsh_element_id;
669 nodetrans[gmsh_node_id],
670 static_cast<boundary_id_type>(entity_to_physical_id[
671 std::make_pair(entity_dim, entity_tag)]));
682 max_elem_dimension_seen=1,
683 min_elem_dimension_seen=3;
686 if (elem_dimensions_seen[i])
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));
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>());
714 for (
const auto & pr : gmsh_physicals)
717 int phys_id = pr.first;
718 unsigned phys_dim = pr.second.first;
719 const std::string & phys_name = pr.second.second;
723 if (phys_dim == max_elem_dimension_seen)
728 else if (phys_dim < max_elem_dimension_seen &&
729 !lower_dimensional_blocks.count(cast_int<boundary_id_type>(phys_id)))
747 typedef std::unordered_multimap<dof_id_type, Elem *> provide_container_t;
748 provide_container_t provide_bcs;
752 if (elem->dim() < max_elem_dimension_seen &&
753 !lower_dimensional_blocks.count(elem->subdomain_id()))
761 for (
auto n : elem->node_index_range())
763 elem->subdomain_id());
768 provide_bcs.insert(std::make_pair(elem->key(), elem));
773 if (elem->dim() == max_elem_dimension_seen)
783 for (
auto sn : elem->side_index_range())
784 for (
const auto & pr :
as_range(provide_bcs.equal_range(elem->key(sn))))
788 std::unique_ptr<Elem> side (elem->build_side_ptr(sn));
791 Elem * lower_dim_elem = pr.second;
794 if (*lower_dim_elem == *side)
800 boundary_id_type bid = cast_int<boundary_id_type>(lower_dim_elem->subdomain_id());
808 if (elem->dim() < max_elem_dimension_seen &&
809 !lower_dimensional_blocks.count(elem->subdomain_id()))
824 libmesh_error_msg(
"Stream is bad! Perhaps the file does not exist?");
836 std::ofstream out_stream (
name.c_str());
839 if (!out_stream.good())
840 libmesh_file_error(
name.c_str());
849 const std::vector<Number> & soln,
850 const std::vector<std::string> & names)
852 LOG_SCOPE(
"write_nodal_data()",
"GmshIO");
870 std::size_t n_boundary_faces = 0;
877 out_stream <<
"$MeshFormat\n";
878 out_stream <<
"2.0 0 " <<
sizeof(
Real) <<
'\n';
879 out_stream <<
"$EndMeshFormat\n";
882 out_stream <<
"$Nodes\n";
890 out_stream <<
"$EndNodes\n";
894 out_stream <<
"$Elements\n";
906 libmesh_assert_less_equal (eletype.
nodes.size(), elem->n_nodes());
909 out_stream << elem->id()+1 <<
" ";
918 << static_cast<unsigned int>(elem->subdomain_id())
920 << elem->processor_id()+1
924 if (eletype.
nodes.size() > 0)
925 for (
auto i : elem->node_index_range())
926 out_stream << elem->node_id(eletype.
nodes[i])+1 <<
" ";
929 for (
auto i : elem->node_index_range())
930 out_stream << elem->node_id(i)+1 <<
" ";
946 if (n_boundary_faces)
952 for (
const auto & t : bc_triples)
956 std::unique_ptr<const Elem> side = elem.
build_side_ptr(std::get<1>(t));
964 libmesh_assert_less_equal (eletype.
nodes.size(), side->n_nodes());
967 out_stream << e_id+1 <<
" ";
983 if (eletype.
nodes.size() > 0)
984 for (
auto i : side->node_index_range())
985 out_stream << side->node_id(eletype.
nodes[i])+1 <<
" ";
989 for (
auto i : side->node_index_range())
990 out_stream << side->node_id(i)+1 <<
" ";
1000 out_stream <<
"$EndElements\n";
1007 const std::vector<Number> * v,
1008 const std::vector<std::string> * solution_names)
1015 std::ofstream out_stream(fname.c_str());
1018 if (!out_stream.good())
1019 libmesh_file_error(fname.c_str());
1028 if ((solution_names !=
nullptr) && (v !=
nullptr))
1030 const unsigned int n_vars =
1031 cast_int<unsigned int>(solution_names->size());
1043 out_stream <<
"$PostFormat\n";
1045 out_stream <<
"1.2 1 " <<
sizeof(
double) <<
"\n";
1047 out_stream <<
"1.2 0 " <<
sizeof(double) <<
"\n";
1048 out_stream <<
"$EndPostFormat\n";
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;
1059 const ElemType elemtype = elem->type();
1109 libmesh_error_msg(
"ERROR: Nonexistent element type " << elem->type());
1115 for (
unsigned int ivar=0; ivar <
n_vars; ivar++)
1117 std::string varname = (*solution_names)[ivar];
1125 out_stream <<
"$View\n" << varname <<
" " << 1 <<
"\n";
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 <<
" "
1153 << nb_text2d_chars <<
" "
1155 << nb_text3d_chars <<
"\n";
1161 std::memcpy(buf, &one,
sizeof(
int));
1162 out_stream.write(buf,
sizeof(
int));
1169 std::memcpy(buf, &one,
sizeof(
double));
1170 out_stream.write(buf,
sizeof(
double));
1173 out_stream <<
"1\n";
1178 const unsigned int nv = elem->n_vertices();
1180 for (
unsigned int d=0; d<3; d++)
1182 for (
unsigned int n=0; n < nv; n++)
1184 const Point & vertex = elem->point(n);
1187 #if defined(LIBMESH_DEFAULT_TRIPLE_PRECISION) || defined(LIBMESH_DEFAULT_QUADRUPLE_PRECISION)
1188 libmesh_warning(
"Gmsh binary writes use only double precision!");
1190 double tmp = double(vertex(d));
1191 std::memcpy(buf, &tmp,
sizeof(
double));
1192 out_stream.write(reinterpret_cast<char *>(buf),
sizeof(
double));
1195 out_stream << vertex(d) <<
" ";
1202 for (
unsigned int i=0; i < nv; i++)
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;
1211 std::memcpy(buf, &tmp,
sizeof(
double));
1212 out_stream.write(reinterpret_cast<char *>(buf),
sizeof(
double));
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;
1226 out_stream <<
"$EndView\n";