Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
grid_in.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
19 #include <deal.II/base/utilities.h>
20 
21 #include <deal.II/grid/grid_in.h>
24 #include <deal.II/grid/tria.h>
25 
26 #include <boost/archive/binary_iarchive.hpp>
27 #include <boost/io/ios_state.hpp>
28 #include <boost/property_tree/ptree.hpp>
29 #include <boost/property_tree/xml_parser.hpp>
30 #include <boost/serialization/serialization.hpp>
31 
32 #include <algorithm>
33 #include <cctype>
34 #include <fstream>
35 #include <functional>
36 #include <map>
37 
38 
39 #ifdef DEAL_II_WITH_NETCDF
40 # include <netcdfcpp.h>
41 #endif
42 
43 #ifdef DEAL_II_WITH_ASSIMP
44 # include <assimp/Importer.hpp> // C++ importer interface
45 # include <assimp/postprocess.h> // Post processing flags
46 # include <assimp/scene.h> // Output data structure
47 #endif
48 
49 
51 
52 
53 namespace
54 {
63  template <int spacedim>
64  void
65  assign_1d_boundary_ids(
66  const std::map<unsigned int, types::boundary_id> &boundary_ids,
68  {
69  if (boundary_ids.size() > 0)
70  for (const auto &cell : triangulation.active_cell_iterators())
71  for (unsigned int f : GeometryInfo<1>::face_indices())
72  if (boundary_ids.find(cell->vertex_index(f)) != boundary_ids.end())
73  {
75  cell->at_boundary(f),
76  ExcMessage(
77  "You are trying to prescribe boundary ids on the face "
78  "of a 1d cell (i.e., on a vertex), but this face is not actually at "
79  "the boundary of the mesh. This is not allowed."));
80  cell->face(f)->set_boundary_id(
81  boundary_ids.find(cell->vertex_index(f))->second);
82  }
83  }
84 
85 
86  template <int dim, int spacedim>
87  void
88  assign_1d_boundary_ids(const std::map<unsigned int, types::boundary_id> &,
90  {
91  // we shouldn't get here since boundary ids are not assigned to
92  // vertices except in 1d
93  Assert(dim != 1, ExcInternalError());
94  }
95 } // namespace
96 
97 template <int dim, int spacedim>
99  : tria(nullptr, typeid(*this).name())
101 {}
102 
103 
104 template <int dim, int spacedim>
105 void
107 {
108  tria = &t;
109 }
110 
111 
112 
113 template <int dim, int spacedim>
114 void
116 {
117  std::string line;
118 
119  // verify that the first, third and fourth lines match
120  // expectations. the second line of the file may essentially be
121  // anything the author of the file chose to identify what's in
122  // there, so we just ensure that we can read it
123  {
124  std::string text[4];
125  text[0] = "# vtk DataFile Version 3.0";
126  text[1] = "****";
127  text[2] = "ASCII";
128  text[3] = "DATASET UNSTRUCTURED_GRID";
129 
130  for (unsigned int i = 0; i < 4; ++i)
131  {
132  getline(in, line);
133  if (i != 1)
134  AssertThrow(
135  line.compare(text[i]) == 0,
136  ExcMessage(
137  std::string(
138  "While reading VTK file, failed to find a header line with text <") +
139  text[i] + ">"));
140  }
141  }
142 
144 
145  std::vector<Point<spacedim>> vertices;
146  std::vector<CellData<dim>> cells;
147  SubCellData subcelldata;
148 
149  std::string keyword;
150 
151  in >> keyword;
152 
154 
155  if (keyword == "POINTS")
156  {
157  unsigned int n_vertices;
158  in >> n_vertices;
159 
160  in >> keyword; // float, double, int, char, etc.
161 
162  for (unsigned int vertex = 0; vertex < n_vertices; ++vertex)
163  {
164  // VTK format always specifies vertex coordinates with 3 components
165  Point<3> x;
166  in >> x(0) >> x(1) >> x(2);
167 
168  vertices.emplace_back();
169  for (unsigned int d = 0; d < spacedim; ++d)
170  vertices.back()(d) = x(d);
171  }
172  }
173 
174  else
175  AssertThrow(false,
176  ExcMessage(
177  "While reading VTK file, failed to find POINTS section"));
178 
179  in >> keyword;
180 
181  unsigned int n_geometric_objects = 0;
182  unsigned int n_ints;
183 
184  if (keyword == "CELLS")
185  {
186  in >> n_geometric_objects;
187  in >> n_ints; // Ignore this, since we don't need it.
188 
189  if (dim == 3)
190  {
191  for (unsigned int count = 0; count < n_geometric_objects; count++)
192  {
193  unsigned int type;
194  in >> type;
195 
196  if (type == 8)
197  {
198  // we assume that the file contains first all cells,
199  // and only then any faces or lines
200  AssertThrow(subcelldata.boundary_quads.size() == 0 &&
201  subcelldata.boundary_lines.size() == 0,
203 
204  cells.emplace_back();
205 
206  for (unsigned int j = 0; j < type; j++) // loop to feed data
207  in >> cells.back().vertices[j];
208 
209  cells.back().material_id = 0;
210  }
211 
212  else if (type == 4)
213  {
214  // we assume that the file contains first all cells,
215  // then all faces, and finally all lines
216  AssertThrow(subcelldata.boundary_lines.size() == 0,
218 
219  subcelldata.boundary_quads.emplace_back();
220 
221  for (unsigned int j = 0; j < type;
222  j++) // loop to feed the data to the boundary
223  in >> subcelldata.boundary_quads.back().vertices[j];
224 
225  subcelldata.boundary_quads.back().material_id = 0;
226  }
227  else if (type == 2)
228  {
229  subcelldata.boundary_lines.emplace_back();
230 
231  for (unsigned int j = 0; j < type;
232  j++) // loop to feed the data to the boundary
233  in >> subcelldata.boundary_lines.back().vertices[j];
234 
235  subcelldata.boundary_lines.back().material_id = 0;
236  }
237 
238  else
239  AssertThrow(
240  false,
241  ExcMessage(
242  "While reading VTK file, unknown file type encountered"));
243  }
244  }
245 
246  else if (dim == 2)
247  {
248  for (unsigned int count = 0; count < n_geometric_objects; count++)
249  {
250  unsigned int type;
251  in >> type;
252 
253  if (type == 4)
254  {
255  // we assume that the file contains first all cells,
256  // and only then any faces
257  AssertThrow(subcelldata.boundary_lines.size() == 0,
259 
260  cells.emplace_back();
261 
262  for (unsigned int j = 0; j < type; j++) // loop to feed data
263  in >> cells.back().vertices[j];
264 
265  cells.back().material_id = 0;
266  }
267 
268  else if (type == 2)
269  {
270  // If this is encountered, the pointer comes out of the loop
271  // and starts processing boundaries.
272  subcelldata.boundary_lines.emplace_back();
273 
274  for (unsigned int j = 0; j < type;
275  j++) // loop to feed the data to the boundary
276  {
277  in >> subcelldata.boundary_lines.back().vertices[j];
278  }
279 
280  subcelldata.boundary_lines.back().material_id = 0;
281  }
282 
283  else
284  AssertThrow(
285  false,
286  ExcMessage(
287  "While reading VTK file, unknown cell type encountered"));
288  }
289  }
290  else if (dim == 1)
291  {
292  for (unsigned int count = 0; count < n_geometric_objects; count++)
293  {
294  unsigned int type;
295  in >> type;
296 
297  AssertThrow(
298  type == 2,
299  ExcMessage(
300  "While reading VTK file, unknown cell type encountered"));
301  cells.emplace_back();
302 
303  for (unsigned int j = 0; j < type; j++) // loop to feed data
304  in >> cells.back().vertices[j];
305 
306  cells.back().material_id = 0;
307  }
308  }
309  else
310  AssertThrow(false,
311  ExcMessage(
312  "While reading VTK file, failed to find CELLS section"));
313 
316 
317  in >> keyword;
318 
319  AssertThrow(
320  keyword == "CELL_TYPES",
321  ExcMessage(std::string(
322  "While reading VTK file, missing CELL_TYPES section. Found <" +
323  keyword + "> instead.")));
324 
325  in >> n_ints;
326  AssertThrow(
327  n_ints == n_geometric_objects,
328  ExcMessage("The VTK reader found a CELL_DATA statement "
329  "that lists a total of " +
330  Utilities::int_to_string(n_ints) +
331  " cell data objects, but this needs to "
332  "equal the number of cells (which is " +
333  Utilities::int_to_string(cells.size()) +
334  ") plus the number of quads (" +
335  Utilities::int_to_string(subcelldata.boundary_quads.size()) +
336  " in 3d or the number of lines (" +
337  Utilities::int_to_string(subcelldata.boundary_lines.size()) +
338  ") in 2d."));
339 
340  int tmp_int;
341  for (unsigned int i = 0; i < n_ints; ++i)
342  in >> tmp_int;
343 
344  // Ignore everything up to CELL_DATA
345  while (in >> keyword)
346  if (keyword == "CELL_DATA")
347  {
348  unsigned int n_ids;
349  in >> n_ids;
350 
351  AssertThrow(n_ids == n_geometric_objects,
352  ExcMessage("The VTK reader found a CELL_DATA statement "
353  "that lists a total of " +
354  Utilities::int_to_string(n_ids) +
355  " cell data objects, but this needs to "
356  "equal the number of cells (which is " +
357  Utilities::int_to_string(cells.size()) +
358  ") plus the number of quads (" +
360  subcelldata.boundary_quads.size()) +
361  " in 3d or the number of lines (" +
363  subcelldata.boundary_lines.size()) +
364  ") in 2d."));
365 
366  const std::vector<std::string> data_sets{"MaterialID",
367  "ManifoldID"};
368 
369  for (unsigned int i = 0; i < data_sets.size(); ++i)
370  {
371  // Ignore everything until we get to a SCALARS data set
372  while (in >> keyword)
373  if (keyword == "SCALARS")
374  {
375  // Now see if we know about this type of data set,
376  // if not, just ignore everything till the next SCALARS
377  // keyword
378  std::string field_name;
379  in >> field_name;
380  if (std::find(data_sets.begin(),
381  data_sets.end(),
382  field_name) == data_sets.end())
383  // The data set here is not one of the ones we know, so
384  // keep ignoring everything until the next SCALARS
385  // keyword.
386  continue;
387 
388  // Now we got somewhere. Proceed from here, assert
389  // that the type of the table is int, and ignore the
390  // rest of the line.
391  // SCALARS MaterialID int 1
392  // (the last number is optional)
393  std::string line;
394  std::getline(in, line);
395  AssertThrow(
396  line.substr(1,
397  std::min(static_cast<std::size_t>(3),
398  line.size() - 1)) == "int",
399  ExcMessage(
400  "While reading VTK file, material- and manifold IDs can only have type 'int'."));
401 
402  in >> keyword;
403  AssertThrow(
404  keyword == "LOOKUP_TABLE",
405  ExcMessage(
406  "While reading VTK file, missing keyword 'LOOKUP_TABLE'."));
407 
408  in >> keyword;
409  AssertThrow(
410  keyword == "default",
411  ExcMessage(
412  "While reading VTK file, missing keyword 'default'."));
413 
414  // read material or manifold ids first for all cells,
415  // then for all faces, and finally for all lines. the
416  // assumption that cells come before all faces and
417  // lines has been verified above via an assertion, so
418  // the order used in the following blocks makes sense
419  for (unsigned int i = 0; i < cells.size(); i++)
420  {
421  int id;
422  in >> id;
423  if (field_name == "MaterialID")
424  cells[i].material_id =
425  static_cast<types::material_id>(id);
426  else if (field_name == "ManifoldID")
427  cells[i].manifold_id =
428  static_cast<types::manifold_id>(id);
429  else
430  Assert(false, ExcInternalError());
431  }
432 
433  if (dim == 3)
434  {
435  for (auto &boundary_quad : subcelldata.boundary_quads)
436  {
437  int id;
438  in >> id;
439  if (field_name == "MaterialID")
440  boundary_quad.material_id =
441  static_cast<types::material_id>(id);
442  else if (field_name == "ManifoldID")
443  boundary_quad.manifold_id =
444  static_cast<types::manifold_id>(id);
445  else
446  Assert(false, ExcInternalError());
447  }
448  for (auto &boundary_line : subcelldata.boundary_lines)
449  {
450  int id;
451  in >> id;
452  if (field_name == "MaterialID")
453  boundary_line.material_id =
454  static_cast<types::material_id>(id);
455  else if (field_name == "ManifoldID")
456  boundary_line.manifold_id =
457  static_cast<types::manifold_id>(id);
458  else
459  Assert(false, ExcInternalError());
460  }
461  }
462  else if (dim == 2)
463  {
464  for (auto &boundary_line : subcelldata.boundary_lines)
465  {
466  int id;
467  in >> id;
468  if (field_name == "MaterialID")
469  boundary_line.material_id =
470  static_cast<types::material_id>(id);
471  else if (field_name == "ManifoldID")
472  boundary_line.manifold_id =
473  static_cast<types::manifold_id>(id);
474  else
475  Assert(false, ExcInternalError());
476  }
477  }
478  }
479  }
480  }
481 
482  Assert(subcelldata.check_consistency(dim), ExcInternalError());
483 
484  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
485 
486  if (dim == spacedim)
488  vertices, cells);
489 
491  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
492 
493  return;
494  }
495  else
496  AssertThrow(false,
497  ExcMessage(
498  "While reading VTK file, failed to find CELLS section"));
499 }
500 
501 
502 
503 template <int dim, int spacedim>
504 void
506 {
507  namespace pt = boost::property_tree;
508  pt::ptree tree;
509  pt::read_xml(in, tree);
510  auto section = tree.get_optional<std::string>("VTKFile.dealiiData");
511 
512  AssertThrow(section,
513  ExcMessage(
514  "While reading a VTU file, failed to find dealiiData section. "
515  "Notice that we can only read grid files in .vtu format that "
516  "were created by the deal.II library, using a call to "
517  "GridOut::write_vtu(), where the flag "
518  "GridOutFlags::Vtu::serialize_triangulation is set to true."));
519 
520  const auto decoded =
521  Utilities::decode_base64({section->begin(), section->end() - 1});
522  const auto string_archive =
523  Utilities::decompress({decoded.begin(), decoded.end()});
524  std::istringstream in_stream(string_archive);
525  boost::archive::binary_iarchive ia(in_stream);
526  tria->load(ia, 0);
527 }
528 
529 
530 template <int dim, int spacedim>
531 void
533 {
534  Assert(tria != nullptr, ExcNoTriangulationSelected());
535  Assert((dim == 2) || (dim == 3), ExcNotImplemented());
536 
537  AssertThrow(in, ExcIO());
538  skip_comment_lines(in, '#'); // skip comments (if any) at beginning of file
539 
540  int tmp;
541 
542  AssertThrow(in, ExcIO());
543  in >> tmp;
544  AssertThrow(in, ExcIO());
545  in >> tmp;
546 
547  // section 2411 describes vertices: see
548  // http://www.sdrl.uc.edu/sdrl/referenceinfo/universalfileformats/file-format-storehouse/universal-dataset-number-2411
549  AssertThrow(tmp == 2411, ExcUnknownSectionType(tmp));
550 
551  std::vector<Point<spacedim>> vertices; // vector of vertex coordinates
552  std::map<int, int>
553  vertex_indices; // # vert in unv (key) ---> # vert in deal.II (value)
554 
555  int no_vertex = 0; // deal.II
556 
557  while (tmp != -1) // we do until reach end of 2411
558  {
559  int no; // unv
560  int dummy;
561  double x[3];
562 
563  AssertThrow(in, ExcIO());
564  in >> no;
565 
566  tmp = no;
567  if (tmp == -1)
568  break;
569 
570  in >> dummy >> dummy >> dummy;
571 
572  AssertThrow(in, ExcIO());
573  in >> x[0] >> x[1] >> x[2];
574 
575  vertices.emplace_back();
576 
577  for (unsigned int d = 0; d < spacedim; d++)
578  vertices.back()(d) = x[d];
579 
580  vertex_indices[no] = no_vertex;
581 
582  no_vertex++;
583  }
584 
585  AssertThrow(in, ExcIO());
586  in >> tmp;
587  AssertThrow(in, ExcIO());
588  in >> tmp;
589 
590  // section 2412 describes elements: see
591  // http://www.sdrl.uc.edu/sdrl/referenceinfo/universalfileformats/file-format-storehouse/universal-dataset-number-2412
592  AssertThrow(tmp == 2412, ExcUnknownSectionType(tmp));
593 
594  std::vector<CellData<dim>> cells; // vector of cells
595  SubCellData subcelldata;
596 
597  std::map<int, int>
598  cell_indices; // # cell in unv (key) ---> # cell in deal.II (value)
599  std::map<int, int>
600  line_indices; // # line in unv (key) ---> # line in deal.II (value)
601  std::map<int, int>
602  quad_indices; // # quad in unv (key) ---> # quad in deal.II (value)
603 
604  int no_cell = 0; // deal.II
605  int no_line = 0; // deal.II
606  int no_quad = 0; // deal.II
607 
608  while (tmp != -1) // we do until reach end of 2412
609  {
610  int no; // unv
611  int type;
612  int dummy;
613 
614  AssertThrow(in, ExcIO());
615  in >> no;
616 
617  tmp = no;
618  if (tmp == -1)
619  break;
620 
621  in >> type >> dummy >> dummy >> dummy >> dummy;
622 
623  AssertThrow((type == 11) || (type == 44) || (type == 94) || (type == 115),
624  ExcUnknownElementType(type));
625 
626  if ((((type == 44) || (type == 94)) && (dim == 2)) ||
627  ((type == 115) && (dim == 3))) // cell
628  {
629  cells.emplace_back();
630 
631  AssertThrow(in, ExcIO());
632  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
633  in >> cells.back().vertices[v];
634 
635  cells.back().material_id = 0;
636 
637  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
638  cells.back().vertices[v] = vertex_indices[cells.back().vertices[v]];
639 
640  cell_indices[no] = no_cell;
641 
642  no_cell++;
643  }
644  else if (((type == 11) && (dim == 2)) ||
645  ((type == 11) && (dim == 3))) // boundary line
646  {
647  AssertThrow(in, ExcIO());
648  in >> dummy >> dummy >> dummy;
649 
650  subcelldata.boundary_lines.emplace_back();
651 
652  AssertThrow(in, ExcIO());
653  for (unsigned int &vertex :
654  subcelldata.boundary_lines.back().vertices)
655  in >> vertex;
656 
657  subcelldata.boundary_lines.back().material_id = 0;
658 
659  for (unsigned int &vertex :
660  subcelldata.boundary_lines.back().vertices)
661  vertex = vertex_indices[vertex];
662 
663  line_indices[no] = no_line;
664 
665  no_line++;
666  }
667  else if (((type == 44) || (type == 94)) && (dim == 3)) // boundary quad
668  {
669  subcelldata.boundary_quads.emplace_back();
670 
671  AssertThrow(in, ExcIO());
672  for (unsigned int &vertex :
673  subcelldata.boundary_quads.back().vertices)
674  in >> vertex;
675 
676  subcelldata.boundary_quads.back().material_id = 0;
677 
678  for (unsigned int &vertex :
679  subcelldata.boundary_quads.back().vertices)
680  vertex = vertex_indices[vertex];
681 
682  quad_indices[no] = no_quad;
683 
684  no_quad++;
685  }
686  else
687  AssertThrow(false,
688  ExcMessage("Unknown element label <" +
690  "> when running in dim=" +
692  }
693 
694  // note that so far all materials and bcs are explicitly set to 0
695  // if we do not need more info on materials and bcs - this is end of file
696  // if we do - section 2467 or 2477 comes
697 
698  in >> tmp; // tmp can be either -1 or end-of-file
699 
700  if (!in.eof())
701  {
702  AssertThrow(in, ExcIO());
703  in >> tmp;
704 
705  // section 2467 (2477) describes (materials - first and bcs - second) or
706  // (bcs - first and materials - second) - sequence depends on which
707  // group is created first: see
708  // http://www.sdrl.uc.edu/sdrl/referenceinfo/universalfileformats/file-format-storehouse/universal-dataset-number-2467
709  AssertThrow((tmp == 2467) || (tmp == 2477), ExcUnknownSectionType(tmp));
710 
711  while (tmp != -1) // we do until reach end of 2467 or 2477
712  {
713  int n_entities; // number of entities in group
714  int id; // id is either material or bc
715  int no; // unv
716  int dummy;
717 
718  AssertThrow(in, ExcIO());
719  in >> dummy;
720 
721  tmp = dummy;
722  if (tmp == -1)
723  break;
724 
725  in >> dummy >> dummy >> dummy >> dummy >> dummy >> dummy >>
726  n_entities;
727 
728  AssertThrow(in, ExcIO());
729  in >> id;
730 
731  const unsigned int n_lines =
732  (n_entities % 2 == 0) ? (n_entities / 2) : ((n_entities + 1) / 2);
733 
734  for (unsigned int line = 0; line < n_lines; line++)
735  {
736  unsigned int n_fragments;
737 
738  if (line == n_lines - 1)
739  n_fragments = (n_entities % 2 == 0) ? (2) : (1);
740  else
741  n_fragments = 2;
742 
743  for (unsigned int no_fragment = 0; no_fragment < n_fragments;
744  no_fragment++)
745  {
746  AssertThrow(in, ExcIO());
747  in >> dummy >> no >> dummy >> dummy;
748 
749  if (cell_indices.count(no) > 0) // cell - material
750  cells[cell_indices[no]].material_id = id;
751 
752  if (line_indices.count(no) > 0) // boundary line - bc
753  subcelldata.boundary_lines[line_indices[no]].material_id =
754  id;
755 
756  if (quad_indices.count(no) > 0) // boundary quad - bc
757  subcelldata.boundary_quads[quad_indices[no]].material_id =
758  id;
759  }
760  }
761  }
762  }
763 
764  Assert(subcelldata.check_consistency(dim), ExcInternalError());
765 
766  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
767 
768  if (dim == spacedim)
770  cells);
771 
773 
774  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
775 }
776 
777 
778 
779 template <int dim, int spacedim>
780 void
782  const bool apply_all_indicators_to_manifolds)
783 {
784  Assert(tria != nullptr, ExcNoTriangulationSelected());
785  AssertThrow(in, ExcIO());
786 
787  // skip comments at start of file
788  skip_comment_lines(in, '#');
789 
790 
791  unsigned int n_vertices;
792  unsigned int n_cells;
793  int dummy;
794 
795  in >> n_vertices >> n_cells >> dummy // number of data vectors
796  >> dummy // cell data
797  >> dummy; // model data
798  AssertThrow(in, ExcIO());
799 
800  // set up array of vertices
801  std::vector<Point<spacedim>> vertices(n_vertices);
802  // set up mapping between numbering
803  // in ucd-file (key) and in the
804  // vertices vector
805  std::map<int, int> vertex_indices;
806 
807  for (unsigned int vertex = 0; vertex < n_vertices; ++vertex)
808  {
809  int vertex_number;
810  double x[3];
811 
812  // read vertex
813  AssertThrow(in, ExcIO());
814  in >> vertex_number >> x[0] >> x[1] >> x[2];
815 
816  // store vertex
817  for (unsigned int d = 0; d < spacedim; ++d)
818  vertices[vertex](d) = x[d];
819  // store mapping; note that
820  // vertices_indices[i] is automatically
821  // created upon first usage
822  vertex_indices[vertex_number] = vertex;
823  }
824 
825  // set up array of cells
826  std::vector<CellData<dim>> cells;
827  SubCellData subcelldata;
828 
829  for (unsigned int cell = 0; cell < n_cells; ++cell)
830  {
831  // note that since in the input
832  // file we found the number of
833  // cells at the top, there
834  // should still be input here,
835  // so check this:
836  AssertThrow(in, ExcIO());
837 
838  std::string cell_type;
839 
840  // we use an unsigned int because we
841  // fill this variable through an read-in process
842  unsigned int material_id;
843 
844  in >> dummy // cell number
845  >> material_id;
846  in >> cell_type;
847 
848  if (((cell_type == "line") && (dim == 1)) ||
849  ((cell_type == "quad") && (dim == 2)) ||
850  ((cell_type == "hex") && (dim == 3)))
851  // found a cell
852  {
853  // allocate and read indices
854  cells.emplace_back();
855  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
856  in >> cells.back().vertices[i];
857 
858  // to make sure that the cast won't fail
861  0,
863  // we use only material_ids in the range from 0 to
864  // numbers::invalid_material_id-1
866 
867  if (apply_all_indicators_to_manifolds)
868  cells.back().manifold_id =
869  static_cast<types::manifold_id>(material_id);
870  cells.back().material_id = material_id;
871 
872  // transform from ucd to
873  // consecutive numbering
874  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
875  if (vertex_indices.find(cells.back().vertices[i]) !=
876  vertex_indices.end())
877  // vertex with this index exists
878  cells.back().vertices[i] =
879  vertex_indices[cells.back().vertices[i]];
880  else
881  {
882  // no such vertex index
883  AssertThrow(false,
885  cells.back().vertices[i]));
886 
887  cells.back().vertices[i] = numbers::invalid_unsigned_int;
888  }
889  }
890  else if ((cell_type == "line") && ((dim == 2) || (dim == 3)))
891  // boundary info
892  {
893  subcelldata.boundary_lines.emplace_back();
894  in >> subcelldata.boundary_lines.back().vertices[0] >>
895  subcelldata.boundary_lines.back().vertices[1];
896 
897  // to make sure that the cast won't fail
900  0,
902  // we use only boundary_ids in the range from 0 to
903  // numbers::internal_face_boundary_id-1
905 
906  // Make sure to set both manifold id and boundary id appropriately in
907  // both cases:
908  // numbers::internal_face_boundary_id and numbers::flat_manifold_id
909  // are ignored in Triangulation::create_triangulation.
910  if (apply_all_indicators_to_manifolds)
911  {
912  subcelldata.boundary_lines.back().boundary_id =
914  subcelldata.boundary_lines.back().manifold_id =
915  static_cast<types::manifold_id>(material_id);
916  }
917  else
918  {
919  subcelldata.boundary_lines.back().boundary_id =
920  static_cast<types::boundary_id>(material_id);
921  subcelldata.boundary_lines.back().manifold_id =
923  }
924 
925  // transform from ucd to
926  // consecutive numbering
927  for (unsigned int &vertex :
928  subcelldata.boundary_lines.back().vertices)
929  if (vertex_indices.find(vertex) != vertex_indices.end())
930  // vertex with this index exists
931  vertex = vertex_indices[vertex];
932  else
933  {
934  // no such vertex index
935  AssertThrow(false, ExcInvalidVertexIndex(cell, vertex));
937  }
938  }
939  else if ((cell_type == "quad") && (dim == 3))
940  // boundary info
941  {
942  subcelldata.boundary_quads.emplace_back();
943  in >> subcelldata.boundary_quads.back().vertices[0] >>
944  subcelldata.boundary_quads.back().vertices[1] >>
945  subcelldata.boundary_quads.back().vertices[2] >>
946  subcelldata.boundary_quads.back().vertices[3];
947 
948  // to make sure that the cast won't fail
951  0,
953  // we use only boundary_ids in the range from 0 to
954  // numbers::internal_face_boundary_id-1
956 
957  // Make sure to set both manifold id and boundary id appropriately in
958  // both cases:
959  // numbers::internal_face_boundary_id and numbers::flat_manifold_id
960  // are ignored in Triangulation::create_triangulation.
961  if (apply_all_indicators_to_manifolds)
962  {
963  subcelldata.boundary_quads.back().boundary_id =
965  subcelldata.boundary_quads.back().manifold_id =
966  static_cast<types::manifold_id>(material_id);
967  }
968  else
969  {
970  subcelldata.boundary_quads.back().boundary_id =
971  static_cast<types::boundary_id>(material_id);
972  subcelldata.boundary_quads.back().manifold_id =
974  }
975 
976  // transform from ucd to
977  // consecutive numbering
978  for (unsigned int &vertex :
979  subcelldata.boundary_quads.back().vertices)
980  if (vertex_indices.find(vertex) != vertex_indices.end())
981  // vertex with this index exists
982  vertex = vertex_indices[vertex];
983  else
984  {
985  // no such vertex index
986  Assert(false, ExcInvalidVertexIndex(cell, vertex));
988  }
989  }
990  else
991  // cannot read this
992  AssertThrow(false, ExcUnknownIdentifier(cell_type));
993  }
994 
995 
996  // check that no forbidden arrays are used
997  Assert(subcelldata.check_consistency(dim), ExcInternalError());
998 
999  AssertThrow(in, ExcIO());
1000 
1001  // do some clean-up on vertices...
1002  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
1003  // ... and cells
1004  if (dim == spacedim)
1006  cells);
1008  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
1009 }
1010 
1011 namespace
1012 {
1013  template <int dim, int spacedim>
1014  class Abaqus_to_UCD
1015  {
1016  public:
1017  Abaqus_to_UCD();
1018 
1019  void
1020  read_in_abaqus(std::istream &in);
1021  void
1022  write_out_avs_ucd(std::ostream &out) const;
1023 
1024  private:
1025  const double tolerance;
1026 
1027  std::vector<double>
1028  get_global_node_numbers(const int face_cell_no,
1029  const int face_cell_face_no) const;
1030 
1031  // NL: Stored as [ global node-id (int), x-coord, y-coord, z-coord ]
1032  std::vector<std::vector<double>> node_list;
1033  // CL: Stored as [ material-id (int), node1, node2, node3, node4, node5,
1034  // node6, node7, node8 ]
1035  std::vector<std::vector<double>> cell_list;
1036  // FL: Stored as [ sideset-id (int), node1, node2, node3, node4 ]
1037  std::vector<std::vector<double>> face_list;
1038  // ELSET: Stored as [ (std::string) elset_name = (std::vector) of cells
1039  // numbers]
1040  std::map<std::string, std::vector<int>> elsets_list;
1041  };
1042 } // namespace
1043 
1044 template <int dim, int spacedim>
1045 void
1047  const bool apply_all_indicators_to_manifolds)
1048 {
1049  Assert(tria != nullptr, ExcNoTriangulationSelected());
1050  // This implementation has only been verified for:
1051  // - 2d grids with codimension 0
1052  // - 3d grids with codimension 0
1053  // - 3d grids with codimension 1
1054  Assert((spacedim == 2 && dim == spacedim) ||
1055  (spacedim == 3 && (dim == spacedim || dim == spacedim - 1)),
1056  ExcNotImplemented());
1057  AssertThrow(in, ExcIO());
1058 
1059  // Read in the Abaqus file into an intermediate object
1060  // that is to be passed along to the UCD reader
1061  Abaqus_to_UCD<dim, spacedim> abaqus_to_ucd;
1062  abaqus_to_ucd.read_in_abaqus(in);
1063 
1064  std::stringstream in_ucd;
1065  abaqus_to_ucd.write_out_avs_ucd(in_ucd);
1066 
1067  // This next call is wrapped in a try-catch for the following reason:
1068  // It ensures that if the Abaqus mesh is read in correctly but produces
1069  // an erroneous result then the user is alerted to the source of the problem
1070  // and doesn't think that they've somehow called the wrong function.
1071  try
1072  {
1073  read_ucd(in_ucd, apply_all_indicators_to_manifolds);
1074  }
1075  catch (std::exception &exc)
1076  {
1077  std::cerr << "Exception on processing internal UCD data: " << std::endl
1078  << exc.what() << std::endl;
1079 
1080  AssertThrow(
1081  false,
1082  ExcMessage(
1083  "Internal conversion from ABAQUS file to UCD format was unsuccessful. "
1084  "More information is provided in an error message printed above. "
1085  "Are you sure that your ABAQUS mesh file conforms with the requirements "
1086  "listed in the documentation?"));
1087  }
1088  catch (...)
1089  {
1090  AssertThrow(
1091  false,
1092  ExcMessage(
1093  "Internal conversion from ABAQUS file to UCD format was unsuccessful. "
1094  "Are you sure that your ABAQUS mesh file conforms with the requirements "
1095  "listed in the documentation?"));
1096  }
1097 }
1098 
1099 
1100 template <int dim, int spacedim>
1101 void
1103 {
1104  Assert(tria != nullptr, ExcNoTriangulationSelected());
1105  Assert(dim == 2, ExcNotImplemented());
1106 
1107  AssertThrow(in, ExcIO());
1108 
1109  // skip comments at start of file
1110  skip_comment_lines(in, '#');
1111 
1112  // first read in identifier string
1113  std::string line;
1114  getline(in, line);
1115 
1116  AssertThrow(line == "MeshVersionFormatted 0", ExcInvalidDBMESHInput(line));
1117 
1118  skip_empty_lines(in);
1119 
1120  // next read dimension
1121  getline(in, line);
1122  AssertThrow(line == "Dimension", ExcInvalidDBMESHInput(line));
1123  unsigned int dimension;
1124  in >> dimension;
1125  AssertThrow(dimension == dim, ExcDBMESHWrongDimension(dimension));
1126  skip_empty_lines(in);
1127 
1128  // now there are a lot of fields of
1129  // which we don't know the exact
1130  // meaning and which are far from
1131  // being properly documented in the
1132  // manual. we skip everything until
1133  // we find a comment line with the
1134  // string "# END". at some point in
1135  // the future, someone may have the
1136  // knowledge to parse and interpret
1137  // the other fields in between as
1138  // well...
1139  while (getline(in, line), line.find("# END") == std::string::npos)
1140  ;
1141  skip_empty_lines(in);
1142 
1143 
1144  // now read vertices
1145  getline(in, line);
1146  AssertThrow(line == "Vertices", ExcInvalidDBMESHInput(line));
1147 
1148  unsigned int n_vertices;
1149  double dummy;
1150 
1151  in >> n_vertices;
1152  std::vector<Point<spacedim>> vertices(n_vertices);
1153  for (unsigned int vertex = 0; vertex < n_vertices; ++vertex)
1154  {
1155  // read vertex coordinates
1156  for (unsigned int d = 0; d < dim; ++d)
1157  in >> vertices[vertex][d];
1158  // read Ref phi_i, whatever that may be
1159  in >> dummy;
1160  }
1161  AssertThrow(in, ExcInvalidDBMeshFormat());
1162 
1163  skip_empty_lines(in);
1164 
1165  // read edges. we ignore them at
1166  // present, so just read them and
1167  // discard the input
1168  getline(in, line);
1169  AssertThrow(line == "Edges", ExcInvalidDBMESHInput(line));
1170 
1171  unsigned int n_edges;
1172  in >> n_edges;
1173  for (unsigned int edge = 0; edge < n_edges; ++edge)
1174  {
1175  // read vertex indices
1176  in >> dummy >> dummy;
1177  // read Ref phi_i, whatever that may be
1178  in >> dummy;
1179  }
1180  AssertThrow(in, ExcInvalidDBMeshFormat());
1181 
1182  skip_empty_lines(in);
1183 
1184 
1185 
1186  // read cracked edges (whatever
1187  // that may be). we ignore them at
1188  // present, so just read them and
1189  // discard the input
1190  getline(in, line);
1191  AssertThrow(line == "CrackedEdges", ExcInvalidDBMESHInput(line));
1192 
1193  in >> n_edges;
1194  for (unsigned int edge = 0; edge < n_edges; ++edge)
1195  {
1196  // read vertex indices
1197  in >> dummy >> dummy;
1198  // read Ref phi_i, whatever that may be
1199  in >> dummy;
1200  }
1201  AssertThrow(in, ExcInvalidDBMeshFormat());
1202 
1203  skip_empty_lines(in);
1204 
1205 
1206  // now read cells.
1207  // set up array of cells
1208  getline(in, line);
1209  AssertThrow(line == "Quadrilaterals", ExcInvalidDBMESHInput(line));
1210 
1211  std::vector<CellData<dim>> cells;
1212  SubCellData subcelldata;
1213  unsigned int n_cells;
1214  in >> n_cells;
1215  for (unsigned int cell = 0; cell < n_cells; ++cell)
1216  {
1217  // read in vertex numbers. they
1218  // are 1-based, so subtract one
1219  cells.emplace_back();
1220  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
1221  {
1222  in >> cells.back().vertices[i];
1223 
1224  AssertThrow((cells.back().vertices[i] >= 1) &&
1225  (static_cast<unsigned int>(cells.back().vertices[i]) <=
1226  vertices.size()),
1227  ExcInvalidVertexIndex(cell, cells.back().vertices[i]));
1228 
1229  --cells.back().vertices[i];
1230  }
1231 
1232  // read and discard Ref phi_i
1233  in >> dummy;
1234  }
1235  AssertThrow(in, ExcInvalidDBMeshFormat());
1236 
1237  skip_empty_lines(in);
1238 
1239 
1240  // then there are again a whole lot
1241  // of fields of which I have no
1242  // clue what they mean. skip them
1243  // all and leave the interpretation
1244  // to other implementors...
1245  while (getline(in, line), ((line.find("End") == std::string::npos) && (in)))
1246  ;
1247  // ok, so we are not at the end of
1248  // the file, that's it, mostly
1249 
1250 
1251  // check that no forbidden arrays are used
1252  Assert(subcelldata.check_consistency(dim), ExcInternalError());
1253 
1254  AssertThrow(in, ExcIO());
1255 
1256  // do some clean-up on vertices...
1257  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
1258  // ...and cells
1260  cells);
1262  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
1263 }
1264 
1265 
1266 
1267 template <int dim, int spacedim>
1268 void
1270 {
1271  Assert(false, ExcNotImplemented());
1272 }
1273 
1274 
1275 
1276 template <>
1277 void
1278 GridIn<2>::read_xda(std::istream &in)
1279 {
1280  Assert(tria != nullptr, ExcNoTriangulationSelected());
1281  AssertThrow(in, ExcIO());
1282 
1283  std::string line;
1284  // skip comments at start of file
1285  getline(in, line);
1286 
1287 
1288  unsigned int n_vertices;
1289  unsigned int n_cells;
1290 
1291  // read cells, throw away rest of line
1292  in >> n_cells;
1293  getline(in, line);
1294 
1295  in >> n_vertices;
1296  getline(in, line);
1297 
1298  // ignore following 8 lines
1299  for (unsigned int i = 0; i < 8; ++i)
1300  getline(in, line);
1301 
1302  // set up array of cells
1303  std::vector<CellData<2>> cells(n_cells);
1304  SubCellData subcelldata;
1305 
1306  for (unsigned int cell = 0; cell < n_cells; ++cell)
1307  {
1308  // note that since in the input
1309  // file we found the number of
1310  // cells at the top, there
1311  // should still be input here,
1312  // so check this:
1313  AssertThrow(in, ExcIO());
1315 
1316  for (unsigned int &vertex : cells[cell].vertices)
1317  in >> vertex;
1318  }
1319 
1320 
1321 
1322  // set up array of vertices
1323  std::vector<Point<2>> vertices(n_vertices);
1324  for (unsigned int vertex = 0; vertex < n_vertices; ++vertex)
1325  {
1326  double x[3];
1327 
1328  // read vertex
1329  in >> x[0] >> x[1] >> x[2];
1330 
1331  // store vertex
1332  for (unsigned int d = 0; d < 2; ++d)
1333  vertices[vertex](d) = x[d];
1334  }
1335  AssertThrow(in, ExcIO());
1336 
1337  // do some clean-up on vertices...
1338  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
1339  // ... and cells
1342  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
1343 }
1344 
1345 
1346 
1347 template <>
1348 void
1349 GridIn<3>::read_xda(std::istream &in)
1350 {
1351  Assert(tria != nullptr, ExcNoTriangulationSelected());
1352  AssertThrow(in, ExcIO());
1353 
1354  static const unsigned int xda_to_dealII_map[] = {0, 1, 5, 4, 3, 2, 6, 7};
1355 
1356  std::string line;
1357  // skip comments at start of file
1358  getline(in, line);
1359 
1360 
1361  unsigned int n_vertices;
1362  unsigned int n_cells;
1363 
1364  // read cells, throw away rest of line
1365  in >> n_cells;
1366  getline(in, line);
1367 
1368  in >> n_vertices;
1369  getline(in, line);
1370 
1371  // ignore following 8 lines
1372  for (unsigned int i = 0; i < 8; ++i)
1373  getline(in, line);
1374 
1375  // set up array of cells
1376  std::vector<CellData<3>> cells(n_cells);
1377  SubCellData subcelldata;
1378 
1379  for (unsigned int cell = 0; cell < n_cells; ++cell)
1380  {
1381  // note that since in the input
1382  // file we found the number of
1383  // cells at the top, there
1384  // should still be input here,
1385  // so check this:
1386  AssertThrow(in, ExcIO());
1388 
1389  unsigned int xda_ordered_nodes[8];
1390 
1391  for (unsigned int &xda_ordered_node : xda_ordered_nodes)
1392  in >> xda_ordered_node;
1393 
1394  for (unsigned int i = 0; i < 8; i++)
1395  cells[cell].vertices[i] = xda_ordered_nodes[xda_to_dealII_map[i]];
1396  }
1397 
1398 
1399 
1400  // set up array of vertices
1401  std::vector<Point<3>> vertices(n_vertices);
1402  for (unsigned int vertex = 0; vertex < n_vertices; ++vertex)
1403  {
1404  double x[3];
1405 
1406  // read vertex
1407  in >> x[0] >> x[1] >> x[2];
1408 
1409  // store vertex
1410  for (unsigned int d = 0; d < 3; ++d)
1411  vertices[vertex](d) = x[d];
1412  }
1413  AssertThrow(in, ExcIO());
1414 
1415  // do some clean-up on vertices...
1416  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
1417  // ... and cells
1420  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
1421 }
1422 
1423 
1424 
1425 template <int dim, int spacedim>
1426 void
1428 {
1429  Assert(tria != nullptr, ExcNoTriangulationSelected());
1430  AssertThrow(in, ExcIO());
1431 
1432  unsigned int n_vertices;
1433  unsigned int n_cells;
1434  unsigned int dummy;
1435  std::string line;
1436  // This array stores maps from the 'entities' to the 'physical tags' for
1437  // points, curves, surfaces and volumes. We use this information later to
1438  // assign boundary ids.
1439  std::array<std::map<int, int>, 4> tag_maps;
1440 
1441  in >> line;
1442 
1443  // first determine file format
1444  unsigned int gmsh_file_format = 0;
1445  if (line == "@f$NOD")
1446  gmsh_file_format = 10;
1447  else if (line == "@f$MeshFormat")
1448  gmsh_file_format = 20;
1449  else
1450  AssertThrow(false, ExcInvalidGMSHInput(line));
1451 
1452  // if file format is 2.0 or greater then we also have to read the rest of the
1453  // header
1454  if (gmsh_file_format == 20)
1455  {
1456  double version;
1457  unsigned int file_type, data_size;
1458 
1459  in >> version >> file_type >> data_size;
1460 
1461  Assert((version >= 2.0) && (version <= 4.1), ExcNotImplemented());
1462  gmsh_file_format = static_cast<unsigned int>(version * 10);
1463 
1464  Assert(file_type == 0, ExcNotImplemented());
1465  Assert(data_size == sizeof(double), ExcNotImplemented());
1466 
1467  // read the end of the header and the first line of the nodes description
1468  // to synch ourselves with the format 1 handling above
1469  in >> line;
1470  AssertThrow(line == "@f$EndMeshFormat", ExcInvalidGMSHInput(line));
1471 
1472  in >> line;
1473  // if the next block is of kind @f$PhysicalNames, ignore it
1474  if (line == "@f$PhysicalNames")
1475  {
1476  do
1477  {
1478  in >> line;
1479  }
1480  while (line != "@f$EndPhysicalNames");
1481  in >> line;
1482  }
1483 
1484  // if the next block is of kind @f$Entities, parse it
1485  if (line == "@f$Entities")
1486  {
1487  unsigned long n_points, n_curves, n_surfaces, n_volumes;
1488 
1489  in >> n_points >> n_curves >> n_surfaces >> n_volumes;
1490  for (unsigned int i = 0; i < n_points; ++i)
1491  {
1492  // parse point ids
1493  int tag;
1494  unsigned int n_physicals;
1495  double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
1496  box_max_z;
1497 
1498  // we only care for 'tag' as key for tag_maps[0]
1499  if (gmsh_file_format > 40)
1500  {
1501  in >> tag >> box_min_x >> box_min_y >> box_min_z >>
1502  n_physicals;
1503  box_max_x = box_min_x;
1504  box_max_y = box_min_y;
1505  box_max_z = box_min_z;
1506  }
1507  else
1508  {
1509  in >> tag >> box_min_x >> box_min_y >> box_min_z >>
1510  box_max_x >> box_max_y >> box_max_z >> n_physicals;
1511  }
1512  // if there is a physical tag, we will use it as boundary id below
1513  AssertThrow(n_physicals < 2,
1514  ExcMessage("More than one tag is not supported!"));
1515  // if there is no physical tag, use 0 as default
1516  int physical_tag = 0;
1517  for (unsigned int j = 0; j < n_physicals; ++j)
1518  in >> physical_tag;
1519  tag_maps[0][tag] = physical_tag;
1520  }
1521  for (unsigned int i = 0; i < n_curves; ++i)
1522  {
1523  // parse curve ids
1524  int tag;
1525  unsigned int n_physicals;
1526  double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
1527  box_max_z;
1528 
1529  // we only care for 'tag' as key for tag_maps[1]
1530  in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
1531  box_max_y >> box_max_z >> n_physicals;
1532  // if there is a physical tag, we will use it as boundary id below
1533  AssertThrow(n_physicals < 2,
1534  ExcMessage("More than one tag is not supported!"));
1535  // if there is no physical tag, use 0 as default
1536  int physical_tag = 0;
1537  for (unsigned int j = 0; j < n_physicals; ++j)
1538  in >> physical_tag;
1539  tag_maps[1][tag] = physical_tag;
1540  // we don't care about the points associated to a curve, but have
1541  // to parse them anyway because their format is unstructured
1542  in >> n_points;
1543  for (unsigned int j = 0; j < n_points; ++j)
1544  in >> tag;
1545  }
1546 
1547  for (unsigned int i = 0; i < n_surfaces; ++i)
1548  {
1549  // parse surface ids
1550  int tag;
1551  unsigned int n_physicals;
1552  double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
1553  box_max_z;
1554 
1555  // we only care for 'tag' as key for tag_maps[2]
1556  in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
1557  box_max_y >> box_max_z >> n_physicals;
1558  // if there is a physical tag, we will use it as boundary id below
1559  AssertThrow(n_physicals < 2,
1560  ExcMessage("More than one tag is not supported!"));
1561  // if there is no physical tag, use 0 as default
1562  int physical_tag = 0;
1563  for (unsigned int j = 0; j < n_physicals; ++j)
1564  in >> physical_tag;
1565  tag_maps[2][tag] = physical_tag;
1566  // we don't care about the curves associated to a surface, but
1567  // have to parse them anyway because their format is unstructured
1568  in >> n_curves;
1569  for (unsigned int j = 0; j < n_curves; ++j)
1570  in >> tag;
1571  }
1572  for (unsigned int i = 0; i < n_volumes; ++i)
1573  {
1574  // parse volume ids
1575  int tag;
1576  unsigned int n_physicals;
1577  double box_min_x, box_min_y, box_min_z, box_max_x, box_max_y,
1578  box_max_z;
1579 
1580  // we only care for 'tag' as key for tag_maps[3]
1581  in >> tag >> box_min_x >> box_min_y >> box_min_z >> box_max_x >>
1582  box_max_y >> box_max_z >> n_physicals;
1583  // if there is a physical tag, we will use it as boundary id below
1584  AssertThrow(n_physicals < 2,
1585  ExcMessage("More than one tag is not supported!"));
1586  // if there is no physical tag, use 0 as default
1587  int physical_tag = 0;
1588  for (unsigned int j = 0; j < n_physicals; ++j)
1589  in >> physical_tag;
1590  tag_maps[3][tag] = physical_tag;
1591  // we don't care about the surfaces associated to a volume, but
1592  // have to parse them anyway because their format is unstructured
1593  in >> n_surfaces;
1594  for (unsigned int j = 0; j < n_surfaces; ++j)
1595  in >> tag;
1596  }
1597  in >> line;
1598  AssertThrow(line == "@f$EndEntities", ExcInvalidGMSHInput(line));
1599  in >> line;
1600  }
1601 
1602  // if the next block is of kind @f$PartitionedEntities, ignore it
1603  if (line == "@f$PartitionedEntities")
1604  {
1605  do
1606  {
1607  in >> line;
1608  }
1609  while (line != "@f$EndPartitionedEntities");
1610  in >> line;
1611  }
1612 
1613  // but the next thing should,
1614  // in any case, be the list of
1615  // nodes:
1616  AssertThrow(line == "@f$Nodes", ExcInvalidGMSHInput(line));
1617  }
1618 
1619  // now read the nodes list
1620  int n_entity_blocks = 1;
1621  if (gmsh_file_format > 40)
1622  {
1623  int min_node_tag;
1624  int max_node_tag;
1625  in >> n_entity_blocks >> n_vertices >> min_node_tag >> max_node_tag;
1626  }
1627  else if (gmsh_file_format == 40)
1628  {
1629  in >> n_entity_blocks >> n_vertices;
1630  }
1631  else
1632  in >> n_vertices;
1633  std::vector<Point<spacedim>> vertices(n_vertices);
1634  // set up mapping between numbering
1635  // in msh-file (nod) and in the
1636  // vertices vector
1637  std::map<int, int> vertex_indices;
1638 
1639  {
1640  unsigned int global_vertex = 0;
1641  for (int entity_block = 0; entity_block < n_entity_blocks; ++entity_block)
1642  {
1643  int parametric;
1644  unsigned long numNodes;
1645 
1646  if (gmsh_file_format < 40)
1647  {
1648  numNodes = n_vertices;
1649  parametric = 0;
1650  }
1651  else
1652  {
1653  // for gmsh_file_format 4.1 the order of tag and dim is reversed,
1654  // but we are ignoring both anyway.
1655  int tagEntity, dimEntity;
1656  in >> tagEntity >> dimEntity >> parametric >> numNodes;
1657  }
1658 
1659  std::vector<int> vertex_numbers;
1660  int vertex_number;
1661  if (gmsh_file_format > 40)
1662  for (unsigned long vertex_per_entity = 0;
1663  vertex_per_entity < numNodes;
1664  ++vertex_per_entity)
1665  {
1666  in >> vertex_number;
1667  vertex_numbers.push_back(vertex_number);
1668  }
1669 
1670  for (unsigned long vertex_per_entity = 0; vertex_per_entity < numNodes;
1671  ++vertex_per_entity, ++global_vertex)
1672  {
1673  int vertex_number;
1674  double x[3];
1675 
1676  // read vertex
1677  if (gmsh_file_format > 40)
1678  {
1679  vertex_number = vertex_numbers[vertex_per_entity];
1680  in >> x[0] >> x[1] >> x[2];
1681  }
1682  else
1683  in >> vertex_number >> x[0] >> x[1] >> x[2];
1684 
1685  for (unsigned int d = 0; d < spacedim; ++d)
1686  vertices[global_vertex](d) = x[d];
1687  // store mapping
1688  vertex_indices[vertex_number] = global_vertex;
1689 
1690  // ignore parametric coordinates
1691  if (parametric != 0)
1692  {
1693  double u = 0.;
1694  double v = 0.;
1695  in >> u >> v;
1696  (void)u;
1697  (void)v;
1698  }
1699  }
1700  }
1701  AssertDimension(global_vertex, n_vertices);
1702  }
1703 
1704  // Assert we reached the end of the block
1705  in >> line;
1706  static const std::string end_nodes_marker[] = {"@f$ENDNOD", "@f$EndNodes"};
1707  AssertThrow(line == end_nodes_marker[gmsh_file_format == 10 ? 0 : 1],
1708  ExcInvalidGMSHInput(line));
1709 
1710  // Now read in next bit
1711  in >> line;
1712  static const std::string begin_elements_marker[] = {"@f$ELM", "@f$Elements"};
1713  AssertThrow(line == begin_elements_marker[gmsh_file_format == 10 ? 0 : 1],
1714  ExcInvalidGMSHInput(line));
1715 
1716  // now read the cell list
1717  if (gmsh_file_format > 40)
1718  {
1719  int min_node_tag;
1720  int max_node_tag;
1721  in >> n_entity_blocks >> n_cells >> min_node_tag >> max_node_tag;
1722  }
1723  else if (gmsh_file_format == 40)
1724  {
1725  in >> n_entity_blocks >> n_cells;
1726  }
1727  else
1728  {
1729  n_entity_blocks = 1;
1730  in >> n_cells;
1731  }
1732 
1733  // set up array of cells and subcells (faces). In 1d, there is currently no
1734  // standard way in deal.II to pass boundary indicators attached to individual
1735  // vertices, so do this by hand via the boundary_ids_1d array
1736  std::vector<CellData<dim>> cells;
1737  SubCellData subcelldata;
1738  std::map<unsigned int, types::boundary_id> boundary_ids_1d;
1739 
1740  {
1741  unsigned int global_cell = 0;
1742  for (int entity_block = 0; entity_block < n_entity_blocks; ++entity_block)
1743  {
1744  unsigned int material_id;
1745  unsigned long numElements;
1746  int cell_type;
1747 
1748  if (gmsh_file_format < 40)
1749  {
1750  material_id = 0;
1751  cell_type = 0;
1752  numElements = n_cells;
1753  }
1754  else if (gmsh_file_format == 40)
1755  {
1756  int tagEntity, dimEntity;
1757  in >> tagEntity >> dimEntity >> cell_type >> numElements;
1758  material_id = tag_maps[dimEntity][tagEntity];
1759  }
1760  else
1761  {
1762  // for gmsh_file_format 4.1 the order of tag and dim is reversed,
1763  int tagEntity, dimEntity;
1764  in >> dimEntity >> tagEntity >> cell_type >> numElements;
1765  material_id = tag_maps[dimEntity][tagEntity];
1766  }
1767 
1768  for (unsigned int cell_per_entity = 0; cell_per_entity < numElements;
1769  ++cell_per_entity, ++global_cell)
1770  {
1771  // note that since in the input
1772  // file we found the number of
1773  // cells at the top, there
1774  // should still be input here,
1775  // so check this:
1776  AssertThrow(in, ExcIO());
1777 
1778  unsigned int nod_num;
1779 
1780  /*
1781  For file format version 1, the format of each cell is as follows:
1782  elm-number elm-type reg-phys reg-elem number-of-nodes
1783  node-number-list
1784 
1785  However, for version 2, the format reads like this:
1786  elm-number elm-type number-of-tags < tag > ... node-number-list
1787 
1788  For version 4, we have:
1789  tag(int) numVert(int) ...
1790 
1791  In the following, we will ignore the element number (we simply
1792  enumerate them in the order in which we read them, and we will
1793  take reg-phys (version 1) or the first tag (version 2, if any tag
1794  is given at all) as material id. For version 4, we already read
1795  the material and the cell type in above.
1796  */
1797 
1798  unsigned int elm_number = 0;
1799  if (gmsh_file_format < 40)
1800  {
1801  in >> elm_number // ELM-NUMBER
1802  >> cell_type; // ELM-TYPE
1803  }
1804 
1805  if (gmsh_file_format < 20)
1806  {
1807  in >> material_id // REG-PHYS
1808  >> dummy // reg_elm
1809  >> nod_num;
1810  }
1811  else if (gmsh_file_format < 40)
1812  {
1813  // read the tags; ignore all but the first one which we will
1814  // interpret as the material_id (for cells) or boundary_id
1815  // (for faces)
1816  unsigned int n_tags;
1817  in >> n_tags;
1818  if (n_tags > 0)
1819  in >> material_id;
1820  else
1821  material_id = 0;
1822 
1823  for (unsigned int i = 1; i < n_tags; ++i)
1824  in >> dummy;
1825 
1827  }
1828  else
1829  {
1830  // ignore tag
1831  int tag;
1832  in >> tag;
1834  }
1835 
1836 
1837  /* `ELM-TYPE'
1838  defines the geometrical type of the N-th element:
1839  `1'
1840  Line (2 nodes, 1 edge).
1841 
1842  `3'
1843  Quadrangle (4 nodes, 4 edges).
1844 
1845  `5'
1846  Hexahedron (8 nodes, 12 edges, 6 faces).
1847 
1848  `15'
1849  Point (1 node).
1850  */
1851 
1852  if (((cell_type == 1) && (dim == 1)) ||
1853  ((cell_type == 3) && (dim == 2)) ||
1854  ((cell_type == 5) && (dim == 3)))
1855  // found a cell
1856  {
1858  ExcMessage(
1859  "Number of nodes does not coincide with the "
1860  "number required for this object"));
1861 
1862  // allocate and read indices
1863  cells.emplace_back();
1864  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
1865  in >> cells.back().vertices[i];
1866 
1867  // to make sure that the cast won't fail
1868  Assert(material_id <=
1870  ExcIndexRange(
1871  material_id,
1872  0,
1874  // we use only material_ids in the range from 0 to
1875  // numbers::invalid_material_id-1
1877 
1878  cells.back().material_id = material_id;
1879 
1880  // transform from ucd to
1881  // consecutive numbering
1882  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
1883  {
1884  AssertThrow(
1885  vertex_indices.find(cells.back().vertices[i]) !=
1886  vertex_indices.end(),
1887  ExcInvalidVertexIndexGmsh(cell_per_entity,
1888  elm_number,
1889  cells.back().vertices[i]));
1890 
1891  // vertex with this index exists
1892  cells.back().vertices[i] =
1893  vertex_indices[cells.back().vertices[i]];
1894  }
1895  }
1896  else if ((cell_type == 1) && ((dim == 2) || (dim == 3)))
1897  // boundary info
1898  {
1899  subcelldata.boundary_lines.emplace_back();
1900  in >> subcelldata.boundary_lines.back().vertices[0] >>
1901  subcelldata.boundary_lines.back().vertices[1];
1902 
1903  // to make sure that the cast won't fail
1904  Assert(material_id <=
1906  ExcIndexRange(
1907  material_id,
1908  0,
1910  // we use only boundary_ids in the range from 0 to
1911  // numbers::internal_face_boundary_id-1
1914 
1915  subcelldata.boundary_lines.back().boundary_id =
1916  static_cast<types::boundary_id>(material_id);
1917 
1918  // transform from ucd to
1919  // consecutive numbering
1920  for (unsigned int &vertex :
1921  subcelldata.boundary_lines.back().vertices)
1922  if (vertex_indices.find(vertex) != vertex_indices.end())
1923  // vertex with this index exists
1924  vertex = vertex_indices[vertex];
1925  else
1926  {
1927  // no such vertex index
1928  AssertThrow(false,
1929  ExcInvalidVertexIndex(cell_per_entity,
1930  vertex));
1932  }
1933  }
1934  else if ((cell_type == 3) && (dim == 3))
1935  // boundary info
1936  {
1937  subcelldata.boundary_quads.emplace_back();
1938  in >> subcelldata.boundary_quads.back().vertices[0] >>
1939  subcelldata.boundary_quads.back().vertices[1] >>
1940  subcelldata.boundary_quads.back().vertices[2] >>
1941  subcelldata.boundary_quads.back().vertices[3];
1942 
1943  // to make sure that the cast won't fail
1944  Assert(material_id <=
1946  ExcIndexRange(
1947  material_id,
1948  0,
1950  // we use only boundary_ids in the range from 0 to
1951  // numbers::internal_face_boundary_id-1
1954 
1955  subcelldata.boundary_quads.back().boundary_id =
1956  static_cast<types::boundary_id>(material_id);
1957 
1958  // transform from gmsh to
1959  // consecutive numbering
1960  for (unsigned int &vertex :
1961  subcelldata.boundary_quads.back().vertices)
1962  if (vertex_indices.find(vertex) != vertex_indices.end())
1963  // vertex with this index exists
1964  vertex = vertex_indices[vertex];
1965  else
1966  {
1967  // no such vertex index
1968  Assert(false,
1969  ExcInvalidVertexIndex(cell_per_entity, vertex));
1971  }
1972  }
1973  else if (cell_type == 15)
1974  {
1975  // read the indices of nodes given
1976  unsigned int node_index = 0;
1977  if (gmsh_file_format < 20)
1978  {
1979  for (unsigned int i = 0; i < nod_num; ++i)
1980  in >> node_index;
1981  }
1982  else
1983  {
1984  in >> node_index;
1985  }
1986 
1987  // we only care about boundary indicators assigned to individual
1988  // vertices in 1d (because otherwise the vertices are not faces)
1989  if (dim == 1)
1990  boundary_ids_1d[vertex_indices[node_index]] = material_id;
1991  }
1992  else
1993  // cannot read this, so throw
1994  // an exception. treat
1995  // triangles and tetrahedra
1996  // specially since this
1997  // deserves a more explicit
1998  // error message
1999  {
2000  AssertThrow(cell_type != 2,
2001  ExcMessage("Found triangles while reading a file "
2002  "in gmsh format. deal.II does not "
2003  "support triangles"));
2004  AssertThrow(cell_type != 11,
2005  ExcMessage("Found tetrahedra while reading a file "
2006  "in gmsh format. deal.II does not "
2007  "support tetrahedra"));
2008 
2009  AssertThrow(false, ExcGmshUnsupportedGeometry(cell_type));
2010  }
2011  }
2012  }
2013  AssertDimension(global_cell, n_cells);
2014  }
2015  // Assert we reached the end of the block
2016  in >> line;
2017  static const std::string end_elements_marker[] = {"@f$ENDELM", "@f$EndElements"};
2018  AssertThrow(line == end_elements_marker[gmsh_file_format == 10 ? 0 : 1],
2019  ExcInvalidGMSHInput(line));
2020 
2021  // check that no forbidden arrays are used
2022  Assert(subcelldata.check_consistency(dim), ExcInternalError());
2023 
2024  AssertThrow(in, ExcIO());
2025 
2026  // check that we actually read some
2027  // cells.
2028  AssertThrow(cells.size() > 0, ExcGmshNoCellInformation());
2029 
2030  // do some clean-up on
2031  // vertices...
2032  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
2033  // ... and cells
2034  if (dim == spacedim)
2036  cells);
2038  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
2039 
2040  // in 1d, we also have to attach boundary ids to vertices, which does not
2041  // currently work through the call above
2042  if (dim == 1)
2043  assign_1d_boundary_ids(boundary_ids_1d, *tria);
2044 }
2045 
2046 
2047 template <>
2048 void
2049 GridIn<1>::read_netcdf(const std::string &)
2050 {
2051  AssertThrow(false, ExcImpossibleInDim(1));
2052 }
2053 
2054 template <>
2055 void
2056 GridIn<1, 2>::read_netcdf(const std::string &)
2057 {
2058  AssertThrow(false, ExcImpossibleInDim(1));
2059 }
2060 
2061 
2062 template <>
2063 void
2064 GridIn<1, 3>::read_netcdf(const std::string &)
2065 {
2066  AssertThrow(false, ExcImpossibleInDim(1));
2067 }
2068 
2069 
2070 template <>
2071 void
2072 GridIn<2, 3>::read_netcdf(const std::string &)
2073 {
2074  Assert(false, ExcNotImplemented());
2075 }
2076 
2077 template <>
2078 void
2079 GridIn<2>::read_netcdf(const std::string &filename)
2080 {
2081 #ifndef DEAL_II_WITH_NETCDF
2082  (void)filename;
2083  AssertThrow(false, ExcNeedsNetCDF());
2084 #else
2085  const unsigned int dim = 2;
2086  const unsigned int spacedim = 2;
2087  Assert(tria != nullptr, ExcNoTriangulationSelected());
2088  // this function assumes the TAU
2089  // grid format.
2090  //
2091  // This format stores 2d grids as
2092  // 3d grids. In particular, a 2d
2093  // grid of n_cells quadrilaterals
2094  // in the y=0 plane is duplicated
2095  // to y=1 to build n_cells
2096  // hexaeders. The surface
2097  // quadrilaterals of this 3d grid
2098  // are marked with boundary
2099  // marker. In the following we read
2100  // in all data required, find the
2101  // boundary marker associated with
2102  // the plane y=0, and extract the
2103  // corresponding 2d data to build a
2104  // Triangulation<2>.
2105 
2106  // In the following, we assume that
2107  // the 2d grid lies in the x-z
2108  // plane (y=0). I.e. we choose:
2109  // point[coord]=0, with coord=1
2110  const unsigned int coord = 1;
2111  // Also x-y-z (0-1-2) point
2112  // coordinates will be transformed
2113  // to x-y (x2d-y2d) coordinates.
2114  // With coord=1 as above, we have
2115  // x-z (0-2) -> (x2d-y2d)
2116  const unsigned int x2d = 0;
2117  const unsigned int y2d = 2;
2118  // For the case, the 2d grid lies
2119  // in x-y or y-z plane instead, the
2120  // following code must be extended
2121  // to find the right value for
2122  // coord, and setting x2d and y2d
2123  // accordingly.
2124 
2125  // First, open the file
2126  NcFile nc(filename.c_str());
2127  AssertThrow(nc.is_valid(), ExcIO());
2128 
2129  // then read n_cells
2130  NcDim *elements_dim = nc.get_dim("no_of_elements");
2131  AssertThrow(elements_dim->is_valid(), ExcIO());
2132  const unsigned int n_cells = elements_dim->size();
2133 
2134  // then we read
2135  // int marker(no_of_markers)
2136  NcDim *marker_dim = nc.get_dim("no_of_markers");
2137  AssertThrow(marker_dim->is_valid(), ExcIO());
2138  const unsigned int n_markers = marker_dim->size();
2139 
2140  NcVar *marker_var = nc.get_var("marker");
2141  AssertThrow(marker_var->is_valid(), ExcIO());
2142  AssertThrow(marker_var->num_dims() == 1, ExcIO());
2143  AssertThrow(static_cast<unsigned int>(marker_var->get_dim(0)->size()) ==
2144  n_markers,
2145  ExcIO());
2146 
2147  std::vector<int> marker(n_markers);
2148  // use &* to convert
2149  // vector<int>::iterator to int *
2150  marker_var->get(&*marker.begin(), n_markers);
2151 
2152  // next we read
2153  // int boundarymarker_of_surfaces(
2154  // no_of_surfaceelements)
2155  NcDim *bquads_dim = nc.get_dim("no_of_surfacequadrilaterals");
2156  AssertThrow(bquads_dim->is_valid(), ExcIO());
2157  const unsigned int n_bquads = bquads_dim->size();
2158 
2159  NcVar *bmarker_var = nc.get_var("boundarymarker_of_surfaces");
2160  AssertThrow(bmarker_var->is_valid(), ExcIO());
2161  AssertThrow(bmarker_var->num_dims() == 1, ExcIO());
2162  AssertThrow(static_cast<unsigned int>(bmarker_var->get_dim(0)->size()) ==
2163  n_bquads,
2164  ExcIO());
2165 
2166  std::vector<int> bmarker(n_bquads);
2167  bmarker_var->get(&*bmarker.begin(), n_bquads);
2168 
2169  // for each marker count the
2170  // number of boundary quads
2171  // which carry this marker
2172  std::map<int, unsigned int> n_bquads_per_bmarker;
2173  for (unsigned int i = 0; i < n_markers; ++i)
2174  {
2175  // the markers should all be
2176  // different
2177  AssertThrow(n_bquads_per_bmarker.find(marker[i]) ==
2178  n_bquads_per_bmarker.end(),
2179  ExcIO());
2180 
2181  n_bquads_per_bmarker[marker[i]] =
2182  count(bmarker.begin(), bmarker.end(), marker[i]);
2183  }
2184  // Note: the n_bquads_per_bmarker
2185  // map could be used to find the
2186  // right coord by finding the
2187  // marker0 such that
2188  // a/ n_bquads_per_bmarker[marker0]==n_cells
2189  // b/ point[coord]==0,
2190  // Condition a/ would hold for at
2191  // least two markers, marker0 and
2192  // marker1, whereas b/ would hold
2193  // for marker0 only. For marker1 we
2194  // then had point[coord]=constant
2195  // with e.g. constant=1 or -1
2196 
2197  // next we read
2198  // int points_of_surfacequadrilaterals(
2199  // no_of_surfacequadrilaterals,
2200  // points_per_surfacequadrilateral)
2201  NcDim *quad_vertices_dim = nc.get_dim("points_per_surfacequadrilateral");
2202  AssertThrow(quad_vertices_dim->is_valid(), ExcIO());
2203  const unsigned int vertices_per_quad = quad_vertices_dim->size();
2205  ExcIO());
2206 
2207  NcVar *vertex_indices_var = nc.get_var("points_of_surfacequadrilaterals");
2208  AssertThrow(vertex_indices_var->is_valid(), ExcIO());
2209  AssertThrow(vertex_indices_var->num_dims() == 2, ExcIO());
2210  AssertThrow(static_cast<unsigned int>(
2211  vertex_indices_var->get_dim(0)->size()) == n_bquads,
2212  ExcIO());
2213  AssertThrow(static_cast<unsigned int>(
2214  vertex_indices_var->get_dim(1)->size()) == vertices_per_quad,
2215  ExcIO());
2216 
2217  std::vector<int> vertex_indices(n_bquads * vertices_per_quad);
2218  vertex_indices_var->get(&*vertex_indices.begin(),
2219  n_bquads,
2220  vertices_per_quad);
2221 
2222  for (const int idx : vertex_indices)
2223  AssertThrow(idx >= 0, ExcIO());
2224 
2225  // next we read
2226  // double points_xc(no_of_points)
2227  // double points_yc(no_of_points)
2228  // double points_zc(no_of_points)
2229  NcDim *vertices_dim = nc.get_dim("no_of_points");
2230  AssertThrow(vertices_dim->is_valid(), ExcIO());
2231  const unsigned int n_vertices = vertices_dim->size();
2232 
2233  NcVar *points_xc = nc.get_var("points_xc");
2234  NcVar *points_yc = nc.get_var("points_yc");
2235  NcVar *points_zc = nc.get_var("points_zc");
2236  AssertThrow(points_xc->is_valid(), ExcIO());
2237  AssertThrow(points_yc->is_valid(), ExcIO());
2238  AssertThrow(points_zc->is_valid(), ExcIO());
2239  AssertThrow(points_xc->num_dims() == 1, ExcIO());
2240  AssertThrow(points_yc->num_dims() == 1, ExcIO());
2241  AssertThrow(points_zc->num_dims() == 1, ExcIO());
2242  AssertThrow(points_yc->get_dim(0)->size() == static_cast<int>(n_vertices),
2243  ExcIO());
2244  AssertThrow(points_zc->get_dim(0)->size() == static_cast<int>(n_vertices),
2245  ExcIO());
2246  AssertThrow(points_xc->get_dim(0)->size() == static_cast<int>(n_vertices),
2247  ExcIO());
2248  std::vector<std::vector<double>> point_values(
2249  3, std::vector<double>(n_vertices));
2250  points_xc->get(point_values[0].data(), n_vertices);
2251  points_yc->get(point_values[1].data(), n_vertices);
2252  points_zc->get(point_values[2].data(), n_vertices);
2253 
2254  // and fill the vertices
2255  std::vector<Point<spacedim>> vertices(n_vertices);
2256  for (unsigned int i = 0; i < n_vertices; ++i)
2257  {
2258  vertices[i](0) = point_values[x2d][i];
2259  vertices[i](1) = point_values[y2d][i];
2260  }
2261 
2262  // For all boundary quads in the
2263  // point[coord]=0 plane add the
2264  // bmarker to zero_plane_markers
2265  std::map<int, bool> zero_plane_markers;
2266  for (unsigned int quad = 0; quad < n_bquads; ++quad)
2267  {
2268  bool zero_plane = true;
2269  for (unsigned int i = 0; i < vertices_per_quad; ++i)
2270  if (point_values[coord][vertex_indices[quad * vertices_per_quad + i]] !=
2271  0)
2272  {
2273  zero_plane = false;
2274  break;
2275  }
2276 
2277  if (zero_plane)
2278  zero_plane_markers[bmarker[quad]] = true;
2279  }
2280  unsigned int sum_of_zero_plane_cells = 0;
2281  for (std::map<int, bool>::const_iterator iter = zero_plane_markers.begin();
2282  iter != zero_plane_markers.end();
2283  ++iter)
2284  sum_of_zero_plane_cells += n_bquads_per_bmarker[iter->first];
2285  AssertThrow(sum_of_zero_plane_cells == n_cells, ExcIO());
2286 
2287  // fill cells with all quads
2288  // associated with
2289  // zero_plane_markers
2290  std::vector<CellData<dim>> cells(n_cells);
2291  for (unsigned int quad = 0, cell = 0; quad < n_bquads; ++quad)
2292  {
2293  bool zero_plane = false;
2294  for (std::map<int, bool>::const_iterator iter =
2295  zero_plane_markers.begin();
2296  iter != zero_plane_markers.end();
2297  ++iter)
2298  if (bmarker[quad] == iter->first)
2299  {
2300  zero_plane = true;
2301  break;
2302  }
2303 
2304  if (zero_plane)
2305  {
2306  for (unsigned int i = 0; i < vertices_per_quad; ++i)
2307  {
2308  Assert(
2309  point_values[coord]
2310  [vertex_indices[quad * vertices_per_quad + i]] == 0,
2311  ExcNotImplemented());
2312  cells[cell].vertices[i] =
2313  vertex_indices[quad * vertices_per_quad + i];
2314  }
2315  ++cell;
2316  }
2317  }
2318 
2319  SubCellData subcelldata;
2320  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
2322  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
2323 #endif
2324 }
2325 
2326 
2327 template <>
2328 void
2329 GridIn<3>::read_netcdf(const std::string &filename)
2330 {
2331 #ifndef DEAL_II_WITH_NETCDF
2332  // do something with the function argument
2333  // to make sure it at least looks used,
2334  // even if it is not
2335  (void)filename;
2336  AssertThrow(false, ExcNeedsNetCDF());
2337 #else
2338  const unsigned int dim = 3;
2339  const unsigned int spacedim = 3;
2340  Assert(tria != nullptr, ExcNoTriangulationSelected());
2341  // this function assumes the TAU
2342  // grid format.
2343 
2344  // First, open the file
2345  NcFile nc(filename.c_str());
2346  AssertThrow(nc.is_valid(), ExcIO());
2347 
2348  // then read n_cells
2349  NcDim *elements_dim = nc.get_dim("no_of_elements");
2350  AssertThrow(elements_dim->is_valid(), ExcIO());
2351  const unsigned int n_cells = elements_dim->size();
2352  // and n_hexes
2353  NcDim *hexes_dim = nc.get_dim("no_of_hexaeders");
2354  AssertThrow(hexes_dim->is_valid(), ExcIO());
2355  const unsigned int n_hexes = hexes_dim->size();
2356  AssertThrow(n_hexes == n_cells,
2357  ExcMessage("deal.II can handle purely hexaedral grids, only."));
2358 
2359  // next we read
2360  // int points_of_hexaeders(
2361  // no_of_hexaeders,
2362  // points_per_hexaeder)
2363  NcDim *hex_vertices_dim = nc.get_dim("points_per_hexaeder");
2364  AssertThrow(hex_vertices_dim->is_valid(), ExcIO());
2365  const unsigned int vertices_per_hex = hex_vertices_dim->size();
2367  ExcIO());
2368 
2369  NcVar *vertex_indices_var = nc.get_var("points_of_hexaeders");
2370  AssertThrow(vertex_indices_var->is_valid(), ExcIO());
2371  AssertThrow(vertex_indices_var->num_dims() == 2, ExcIO());
2372  AssertThrow(static_cast<unsigned int>(
2373  vertex_indices_var->get_dim(0)->size()) == n_cells,
2374  ExcIO());
2375  AssertThrow(static_cast<unsigned int>(
2376  vertex_indices_var->get_dim(1)->size()) == vertices_per_hex,
2377  ExcIO());
2378 
2379  std::vector<int> vertex_indices(n_cells * vertices_per_hex);
2380  // use &* to convert
2381  // vector<int>::iterator to int *
2382  vertex_indices_var->get(&*vertex_indices.begin(), n_cells, vertices_per_hex);
2383 
2384  for (const int idx : vertex_indices)
2385  AssertThrow(idx >= 0, ExcIO());
2386 
2387  // next we read
2388  // double points_xc(no_of_points)
2389  // double points_yc(no_of_points)
2390  // double points_zc(no_of_points)
2391  NcDim *vertices_dim = nc.get_dim("no_of_points");
2392  AssertThrow(vertices_dim->is_valid(), ExcIO());
2393  const unsigned int n_vertices = vertices_dim->size();
2394 
2395  NcVar *points_xc = nc.get_var("points_xc");
2396  NcVar *points_yc = nc.get_var("points_yc");
2397  NcVar *points_zc = nc.get_var("points_zc");
2398  AssertThrow(points_xc->is_valid(), ExcIO());
2399  AssertThrow(points_yc->is_valid(), ExcIO());
2400  AssertThrow(points_zc->is_valid(), ExcIO());
2401  AssertThrow(points_xc->num_dims() == 1, ExcIO());
2402  AssertThrow(points_yc->num_dims() == 1, ExcIO());
2403  AssertThrow(points_zc->num_dims() == 1, ExcIO());
2404  AssertThrow(points_yc->get_dim(0)->size() == static_cast<int>(n_vertices),
2405  ExcIO());
2406  AssertThrow(points_zc->get_dim(0)->size() == static_cast<int>(n_vertices),
2407  ExcIO());
2408  AssertThrow(points_xc->get_dim(0)->size() == static_cast<int>(n_vertices),
2409  ExcIO());
2410  std::vector<std::vector<double>> point_values(
2411  3, std::vector<double>(n_vertices));
2412  points_xc->get(point_values[0].data(), n_vertices);
2413  points_yc->get(point_values[1].data(), n_vertices);
2414  points_zc->get(point_values[2].data(), n_vertices);
2415 
2416  // and fill the vertices
2417  std::vector<Point<spacedim>> vertices(n_vertices);
2418  for (unsigned int i = 0; i < n_vertices; ++i)
2419  {
2420  vertices[i](0) = point_values[0][i];
2421  vertices[i](1) = point_values[1][i];
2422  vertices[i](2) = point_values[2][i];
2423  }
2424 
2425  // and cells
2426  std::vector<CellData<dim>> cells(n_cells);
2427  for (unsigned int cell = 0; cell < n_cells; ++cell)
2428  for (unsigned int i = 0; i < vertices_per_hex; ++i)
2429  cells[cell].vertices[i] = vertex_indices[cell * vertices_per_hex + i];
2430 
2431  // for setting up the SubCellData
2432  // we read the vertex indices of
2433  // the boundary quadrilaterals and
2434  // their boundary markers
2435 
2436  // first we read
2437  // int points_of_surfacequadrilaterals(
2438  // no_of_surfacequadrilaterals,
2439  // points_per_surfacequadrilateral)
2440  NcDim *quad_vertices_dim = nc.get_dim("points_per_surfacequadrilateral");
2441  AssertThrow(quad_vertices_dim->is_valid(), ExcIO());
2442  const unsigned int vertices_per_quad = quad_vertices_dim->size();
2444  ExcIO());
2445 
2446  NcVar *bvertex_indices_var = nc.get_var("points_of_surfacequadrilaterals");
2447  AssertThrow(bvertex_indices_var->is_valid(), ExcIO());
2448  AssertThrow(bvertex_indices_var->num_dims() == 2, ExcIO());
2449  const unsigned int n_bquads = bvertex_indices_var->get_dim(0)->size();
2450  AssertThrow(static_cast<unsigned int>(
2451  bvertex_indices_var->get_dim(1)->size()) ==
2453  ExcIO());
2454 
2455  std::vector<int> bvertex_indices(n_bquads * vertices_per_quad);
2456  bvertex_indices_var->get(&*bvertex_indices.begin(),
2457  n_bquads,
2458  vertices_per_quad);
2459 
2460  // next we read
2461  // int boundarymarker_of_surfaces(
2462  // no_of_surfaceelements)
2463  NcDim *bquads_dim = nc.get_dim("no_of_surfacequadrilaterals");
2464  AssertThrow(bquads_dim->is_valid(), ExcIO());
2465  AssertThrow(static_cast<unsigned int>(bquads_dim->size()) == n_bquads,
2466  ExcIO());
2467 
2468  NcVar *bmarker_var = nc.get_var("boundarymarker_of_surfaces");
2469  AssertThrow(bmarker_var->is_valid(), ExcIO());
2470  AssertThrow(bmarker_var->num_dims() == 1, ExcIO());
2471  AssertThrow(static_cast<unsigned int>(bmarker_var->get_dim(0)->size()) ==
2472  n_bquads,
2473  ExcIO());
2474 
2475  std::vector<int> bmarker(n_bquads);
2476  bmarker_var->get(&*bmarker.begin(), n_bquads);
2477  // we only handle boundary
2478  // indicators that fit into an
2479  // types::boundary_id. Also, we don't
2480  // take numbers::internal_face_boundary_id
2481  // as it denotes an internal face
2482  for (const int id : bmarker)
2483  {
2484  Assert(0 <= id && static_cast<types::boundary_id>(id) !=
2486  ExcIO());
2487  (void)id;
2488  }
2489 
2490  // finally we setup the boundary
2491  // information
2492  SubCellData subcelldata;
2493  subcelldata.boundary_quads.resize(n_bquads);
2494  for (unsigned int i = 0; i < n_bquads; ++i)
2495  {
2496  for (unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face; ++v)
2497  subcelldata.boundary_quads[i].vertices[v] =
2498  bvertex_indices[i * GeometryInfo<dim>::vertices_per_face + v];
2499  subcelldata.boundary_quads[i].boundary_id =
2500  static_cast<types::boundary_id>(bmarker[i]);
2501  }
2502 
2503  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
2505  cells);
2507  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
2508 #endif
2509 }
2510 
2511 
2512 template <int dim, int spacedim>
2513 void
2515  std::string & header,
2516  std::vector<unsigned int> &tecplot2deal,
2517  unsigned int & n_vars,
2518  unsigned int & n_vertices,
2519  unsigned int & n_cells,
2520  std::vector<unsigned int> &IJK,
2521  bool & structured,
2522  bool & blocked)
2523 {
2524  Assert(tecplot2deal.size() == dim, ExcInternalError());
2525  Assert(IJK.size() == dim, ExcInternalError());
2526  // initialize the output variables
2527  n_vars = 0;
2528  n_vertices = 0;
2529  n_cells = 0;
2530  switch (dim)
2531  {
2532  case 3:
2533  IJK[2] = 0;
2535  case 2:
2536  IJK[1] = 0;
2538  case 1:
2539  IJK[0] = 0;
2540  }
2541  structured = true;
2542  blocked = false;
2543 
2544  // convert the string to upper case
2545  std::transform(header.begin(), header.end(), header.begin(), ::toupper);
2546 
2547  // replace all tabs, commas, newlines by
2548  // whitespaces
2549  std::replace(header.begin(), header.end(), '\t', ' ');
2550  std::replace(header.begin(), header.end(), ',', ' ');
2551  std::replace(header.begin(), header.end(), '\n', ' ');
2552 
2553  // now remove whitespace in front of and
2554  // after '='
2555  std::string::size_type pos = header.find('=');
2556 
2557  while (pos != static_cast<std::string::size_type>(std::string::npos))
2558  if (header[pos + 1] == ' ')
2559  header.erase(pos + 1, 1);
2560  else if (header[pos - 1] == ' ')
2561  {
2562  header.erase(pos - 1, 1);
2563  --pos;
2564  }
2565  else
2566  pos = header.find('=', ++pos);
2567 
2568  // split the string into individual entries
2569  std::vector<std::string> entries =
2570  Utilities::break_text_into_lines(header, 1, ' ');
2571 
2572  // now go through the list and try to extract
2573  for (unsigned int i = 0; i < entries.size(); ++i)
2574  {
2575  if (Utilities::match_at_string_start(entries[i], "VARIABLES=\""))
2576  {
2577  ++n_vars;
2578  // we assume, that the first variable
2579  // is x or no coordinate at all (not y or z)
2580  if (Utilities::match_at_string_start(entries[i], "VARIABLES=\"X\""))
2581  {
2582  tecplot2deal[0] = 0;
2583  }
2584  ++i;
2585  while (entries[i][0] == '"')
2586  {
2587  if (entries[i] == "\"X\"")
2588  tecplot2deal[0] = n_vars;
2589  else if (entries[i] == "\"Y\"")
2590  {
2591  // we assume, that y contains
2592  // zero data in 1d, so do
2593  // nothing
2594  if (dim > 1)
2595  tecplot2deal[1] = n_vars;
2596  }
2597  else if (entries[i] == "\"Z\"")
2598  {
2599  // we assume, that z contains
2600  // zero data in 1d and 2d, so
2601  // do nothing
2602  if (dim > 2)
2603  tecplot2deal[2] = n_vars;
2604  }
2605  ++n_vars;
2606  ++i;
2607  }
2608  // set i back, so that the next
2609  // string is treated correctly
2610  --i;
2611 
2612  AssertThrow(
2613  n_vars >= dim,
2614  ExcMessage(
2615  "Tecplot file must contain at least one variable for each dimension"));
2616  for (unsigned int d = 1; d < dim; ++d)
2617  AssertThrow(
2618  tecplot2deal[d] > 0,
2619  ExcMessage(
2620  "Tecplot file must contain at least one variable for each dimension."));
2621  }
2622  else if (Utilities::match_at_string_start(entries[i], "ZONETYPE=ORDERED"))
2623  structured = true;
2624  else if (Utilities::match_at_string_start(entries[i],
2625  "ZONETYPE=FELINESEG") &&
2626  dim == 1)
2627  structured = false;
2628  else if (Utilities::match_at_string_start(entries[i],
2629  "ZONETYPE=FEQUADRILATERAL") &&
2630  dim == 2)
2631  structured = false;
2632  else if (Utilities::match_at_string_start(entries[i],
2633  "ZONETYPE=FEBRICK") &&
2634  dim == 3)
2635  structured = false;
2636  else if (Utilities::match_at_string_start(entries[i], "ZONETYPE="))
2637  // unsupported ZONETYPE
2638  {
2639  AssertThrow(false,
2640  ExcMessage(
2641  "The tecplot file contains an unsupported ZONETYPE."));
2642  }
2643  else if (Utilities::match_at_string_start(entries[i],
2644  "DATAPACKING=POINT"))
2645  blocked = false;
2646  else if (Utilities::match_at_string_start(entries[i],
2647  "DATAPACKING=BLOCK"))
2648  blocked = true;
2649  else if (Utilities::match_at_string_start(entries[i], "F=POINT"))
2650  {
2651  structured = true;
2652  blocked = false;
2653  }
2654  else if (Utilities::match_at_string_start(entries[i], "F=BLOCK"))
2655  {
2656  structured = true;
2657  blocked = true;
2658  }
2659  else if (Utilities::match_at_string_start(entries[i], "F=FEPOINT"))
2660  {
2661  structured = false;
2662  blocked = false;
2663  }
2664  else if (Utilities::match_at_string_start(entries[i], "F=FEBLOCK"))
2665  {
2666  structured = false;
2667  blocked = true;
2668  }
2669  else if (Utilities::match_at_string_start(entries[i],
2670  "ET=QUADRILATERAL") &&
2671  dim == 2)
2672  structured = false;
2673  else if (Utilities::match_at_string_start(entries[i], "ET=BRICK") &&
2674  dim == 3)
2675  structured = false;
2676  else if (Utilities::match_at_string_start(entries[i], "ET="))
2677  // unsupported ElementType
2678  {
2679  AssertThrow(
2680  false,
2681  ExcMessage(
2682  "The tecplot file contains an unsupported ElementType."));
2683  }
2684  else if (Utilities::match_at_string_start(entries[i], "I="))
2685  IJK[0] = Utilities::get_integer_at_position(entries[i], 2).first;
2686  else if (Utilities::match_at_string_start(entries[i], "J="))
2687  {
2688  IJK[1] = Utilities::get_integer_at_position(entries[i], 2).first;
2689  AssertThrow(
2690  dim > 1 || IJK[1] == 1,
2691  ExcMessage(
2692  "Parameter 'J=' found in tecplot, although this is only possible for dimensions greater than 1."));
2693  }
2694  else if (Utilities::match_at_string_start(entries[i], "K="))
2695  {
2696  IJK[2] = Utilities::get_integer_at_position(entries[i], 2).first;
2697  AssertThrow(
2698  dim > 2 || IJK[2] == 1,
2699  ExcMessage(
2700  "Parameter 'K=' found in tecplot, although this is only possible for dimensions greater than 2."));
2701  }
2702  else if (Utilities::match_at_string_start(entries[i], "N="))
2703  n_vertices = Utilities::get_integer_at_position(entries[i], 2).first;
2704  else if (Utilities::match_at_string_start(entries[i], "E="))
2705  n_cells = Utilities::get_integer_at_position(entries[i], 2).first;
2706  }
2707 
2708  // now we have read all the fields we are
2709  // interested in. do some checks and
2710  // calculate the variables
2711  if (structured)
2712  {
2713  n_vertices = 1;
2714  n_cells = 1;
2715  for (unsigned int d = 0; d < dim; ++d)
2716  {
2717  AssertThrow(
2718  IJK[d] > 0,
2719  ExcMessage(
2720  "Tecplot file does not contain a complete and consistent set of parameters"));
2721  n_vertices *= IJK[d];
2722  n_cells *= (IJK[d] - 1);
2723  }
2724  }
2725  else
2726  {
2727  AssertThrow(
2728  n_vertices > 0,
2729  ExcMessage(
2730  "Tecplot file does not contain a complete and consistent set of parameters"));
2731  if (n_cells == 0)
2732  // this means an error, although
2733  // tecplot itself accepts entries like
2734  // 'J=20' instead of 'E=20'. therefore,
2735  // take the max of IJK
2736  n_cells = *std::max_element(IJK.begin(), IJK.end());
2737  AssertThrow(
2738  n_cells > 0,
2739  ExcMessage(
2740  "Tecplot file does not contain a complete and consistent set of parameters"));
2741  }
2742 }
2743 
2744 
2745 
2746 template <>
2747 void
2748 GridIn<2>::read_tecplot(std::istream &in)
2749 {
2750  const unsigned int dim = 2;
2751  const unsigned int spacedim = 2;
2752  Assert(tria != nullptr, ExcNoTriangulationSelected());
2753  AssertThrow(in, ExcIO());
2754 
2755  // skip comments at start of file
2756  skip_comment_lines(in, '#');
2757 
2758  // some strings for parsing the header
2759  std::string line, header;
2760 
2761  // first, concatenate all header lines
2762  // create a searchstring with almost all
2763  // letters. exclude e and E from the letters
2764  // to search, as they might appear in
2765  // exponential notation
2766  std::string letters = "abcdfghijklmnopqrstuvwxyzABCDFGHIJKLMNOPQRSTUVWXYZ";
2767 
2768  getline(in, line);
2769  while (line.find_first_of(letters) != std::string::npos)
2770  {
2771  header += " " + line;
2772  getline(in, line);
2773  }
2774 
2775  // now create some variables holding
2776  // important information on the mesh, get
2777  // this information from the header string
2778  std::vector<unsigned int> tecplot2deal(dim);
2779  std::vector<unsigned int> IJK(dim);
2780  unsigned int n_vars, n_vertices, n_cells;
2781  bool structured, blocked;
2782 
2783  parse_tecplot_header(header,
2784  tecplot2deal,
2785  n_vars,
2786  n_vertices,
2787  n_cells,
2788  IJK,
2789  structured,
2790  blocked);
2791 
2792  // reserve space for vertices. note, that in
2793  // tecplot vertices are ordered beginning
2794  // with 1, whereas in deal all indices start
2795  // with 0. in order not to use -1 for all the
2796  // connectivity information, a 0th vertex
2797  // (unused) is inserted at the origin.
2798  std::vector<Point<spacedim>> vertices(n_vertices + 1);
2799  vertices[0] = Point<spacedim>();
2800  // reserve space for cells
2801  std::vector<CellData<dim>> cells(n_cells);
2802  SubCellData subcelldata;
2803 
2804  if (blocked)
2805  {
2806  // blocked data format. first we get all
2807  // the values of the first variable for
2808  // all points, after that we get all
2809  // values for the second variable and so
2810  // on.
2811 
2812  // dummy variable to read in all the info
2813  // we do not want to use
2814  double dummy;
2815  // which is the first index to read in
2816  // the loop (see below)
2817  unsigned int next_index = 0;
2818 
2819  // note, that we have already read the
2820  // first line containing the first variable
2821  if (tecplot2deal[0] == 0)
2822  {
2823  // we need the information in this
2824  // line, so extract it
2825  std::vector<std::string> first_var =
2827  char *endptr;
2828  for (unsigned int i = 1; i < first_var.size() + 1; ++i)
2829  vertices[i](0) = std::strtod(first_var[i - 1].c_str(), &endptr);
2830 
2831  // if there are many points, the data
2832  // for this var might continue in the
2833  // next line(s)
2834  for (unsigned int j = first_var.size() + 1; j < n_vertices + 1; ++j)
2835  in >> vertices[j](next_index);
2836  // now we got all values of the first
2837  // variable, so increase the counter
2838  next_index = 1;
2839  }
2840 
2841  // main loop over all variables
2842  for (unsigned int i = 1; i < n_vars; ++i)
2843  {
2844  // if we read all the important
2845  // variables and do not want to
2846  // read further, because we are
2847  // using a structured grid, we can
2848  // stop here (and skip, for
2849  // example, a whole lot of solution
2850  // variables)
2851  if (next_index == dim && structured)
2852  break;
2853 
2854  if ((next_index < dim) && (i == tecplot2deal[next_index]))
2855  {
2856  // we need this line, read it in
2857  for (unsigned int j = 1; j < n_vertices + 1; ++j)
2858  in >> vertices[j](next_index);
2859  ++next_index;
2860  }
2861  else
2862  {
2863  // we do not need this line, read
2864  // it in and discard it
2865  for (unsigned int j = 1; j < n_vertices + 1; ++j)
2866  in >> dummy;
2867  }
2868  }
2869  Assert(next_index == dim, ExcInternalError());
2870  }
2871  else
2872  {
2873  // the data is not blocked, so we get all
2874  // the variables for one point, then the
2875  // next and so on. create a vector to
2876  // hold these components
2877  std::vector<double> vars(n_vars);
2878 
2879  // now fill the first vertex. note, that we
2880  // have already read the first line
2881  // containing the first vertex
2882  std::vector<std::string> first_vertex =
2884  char *endptr;
2885  for (unsigned int d = 0; d < dim; ++d)
2886  vertices[1](d) =
2887  std::strtod(first_vertex[tecplot2deal[d]].c_str(), &endptr);
2888 
2889  // read the remaining vertices from the
2890  // list
2891  for (unsigned int v = 2; v < n_vertices + 1; ++v)
2892  {
2893  for (unsigned int i = 0; i < n_vars; ++i)
2894  in >> vars[i];
2895  // fill the vertex
2896  // coordinates. respect the position
2897  // of coordinates in the list of
2898  // variables
2899  for (unsigned int i = 0; i < dim; ++i)
2900  vertices[v](i) = vars[tecplot2deal[i]];
2901  }
2902  }
2903 
2904  if (structured)
2905  {
2906  // this is the part of the code that only
2907  // works in 2d
2908  unsigned int I = IJK[0], J = IJK[1];
2909 
2910  unsigned int cell = 0;
2911  // set up array of cells
2912  for (unsigned int j = 0; j < J - 1; ++j)
2913  for (unsigned int i = 1; i < I; ++i)
2914  {
2915  cells[cell].vertices[0] = i + j * I;
2916  cells[cell].vertices[1] = i + 1 + j * I;
2917  cells[cell].vertices[2] = i + 1 + (j + 1) * I;
2918  cells[cell].vertices[3] = i + (j + 1) * I;
2919  ++cell;
2920  }
2921  Assert(cell == n_cells, ExcInternalError());
2922  std::vector<unsigned int> boundary_vertices(2 * I + 2 * J - 4);
2923  unsigned int k = 0;
2924  for (unsigned int i = 1; i < I + 1; ++i)
2925  {
2926  boundary_vertices[k] = i;
2927  ++k;
2928  boundary_vertices[k] = i + (J - 1) * I;
2929  ++k;
2930  }
2931  for (unsigned int j = 1; j < J - 1; ++j)
2932  {
2933  boundary_vertices[k] = 1 + j * I;
2934  ++k;
2935  boundary_vertices[k] = I + j * I;
2936  ++k;
2937  }
2938  Assert(k == boundary_vertices.size(), ExcInternalError());
2939  // delete the duplicated vertices at the
2940  // boundary, which occur, e.g. in c-type
2941  // or o-type grids around a body
2942  // (airfoil). this automatically deletes
2943  // unused vertices as well.
2945  cells,
2946  subcelldata,
2947  boundary_vertices);
2948  }
2949  else
2950  {
2951  // set up array of cells, unstructured
2952  // mode, so the connectivity is
2953  // explicitly given
2954  for (unsigned int i = 0; i < n_cells; ++i)
2955  {
2956  // note that since in the input file
2957  // we found the number of cells at
2958  // the top, there should still be
2959  // input here, so check this:
2960  AssertThrow(in, ExcIO());
2961 
2962  // get the connectivity from the
2963  // input file. the vertices are
2964  // ordered like in the ucd format
2965  for (unsigned int &vertex : cells[i].vertices)
2966  in >> vertex;
2967  }
2968  // do some clean-up on vertices
2969  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
2970  }
2971 
2972  // check that no forbidden arrays are
2973  // used. as we do not read in any
2974  // subcelldata, nothing should happen here.
2975  Assert(subcelldata.check_consistency(dim), ExcInternalError());
2976  AssertThrow(in, ExcIO());
2977 
2978  // do some cleanup on cells
2980  cells);
2982  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
2983 }
2984 
2985 
2986 
2987 template <int dim, int spacedim>
2988 void
2990 {
2991  Assert(false, ExcNotImplemented());
2992 }
2993 
2994 
2995 
2996 template <int dim, int spacedim>
2997 void
2998 GridIn<dim, spacedim>::read_assimp(const std::string &filename,
2999  const unsigned int mesh_index,
3000  const bool remove_duplicates,
3001  const double tol,
3002  const bool ignore_unsupported_types)
3003 {
3004 #ifdef DEAL_II_WITH_ASSIMP
3005  // Only good for surface grids.
3006  AssertThrow(dim < 3, ExcImpossibleInDim(dim));
3007 
3008  // Create an instance of the Importer class
3009  Assimp::Importer importer;
3010 
3011  // And have it read the given file with some postprocessing
3012  const aiScene *scene =
3013  importer.ReadFile(filename.c_str(),
3014  aiProcess_RemoveComponent |
3015  aiProcess_JoinIdenticalVertices |
3016  aiProcess_ImproveCacheLocality | aiProcess_SortByPType |
3017  aiProcess_OptimizeGraph | aiProcess_OptimizeMeshes);
3018 
3019  // If the import failed, report it
3020  AssertThrow(scene != nullptr, ExcMessage(importer.GetErrorString()));
3021 
3022  AssertThrow(scene->mNumMeshes != 0,
3023  ExcMessage("Input file contains no meshes."));
3024 
3025  AssertThrow((mesh_index == numbers::invalid_unsigned_int) ||
3026  (mesh_index < scene->mNumMeshes),
3027  ExcMessage("Too few meshes in the file."));
3028 
3029  unsigned int start_mesh =
3030  (mesh_index == numbers::invalid_unsigned_int ? 0 : mesh_index);
3031  unsigned int end_mesh =
3032  (mesh_index == numbers::invalid_unsigned_int ? scene->mNumMeshes :
3033  mesh_index + 1);
3034 
3035  // Deal.II objects are created empty, and then filled with imported file.
3036  std::vector<Point<spacedim>> vertices;
3037  std::vector<CellData<dim>> cells;
3038  SubCellData subcelldata;
3039 
3040  // A series of counters to merge cells.
3041  unsigned int v_offset = 0;
3042  unsigned int c_offset = 0;
3043 
3044  // The index of the mesh will be used as a material index.
3045  for (unsigned int m = start_mesh; m < end_mesh; ++m)
3046  {
3047  const aiMesh *mesh = scene->mMeshes[m];
3048 
3049  // Check that we know what to do with this mesh, otherwise just
3050  // ignore it
3051  if ((dim == 2) && mesh->mPrimitiveTypes != aiPrimitiveType_POLYGON)
3052  {
3053  AssertThrow(ignore_unsupported_types,
3054  ExcMessage("Incompatible mesh " + std::to_string(m) +
3055  "/" + std::to_string(scene->mNumMeshes)));
3056  continue;
3057  }
3058  else if ((dim == 1) && mesh->mPrimitiveTypes != aiPrimitiveType_LINE)
3059  {
3060  AssertThrow(ignore_unsupported_types,
3061  ExcMessage("Incompatible mesh " + std::to_string(m) +
3062  "/" + std::to_string(scene->mNumMeshes)));
3063  continue;
3064  }
3065  // Vertices
3066  const unsigned int n_vertices = mesh->mNumVertices;
3067  const aiVector3D * mVertices = mesh->mVertices;
3068 
3069  // Faces
3070  const unsigned int n_faces = mesh->mNumFaces;
3071  const aiFace * mFaces = mesh->mFaces;
3072 
3073  vertices.resize(v_offset + n_vertices);
3074  cells.resize(c_offset + n_faces);
3075 
3076  for (unsigned int i = 0; i < n_vertices; ++i)
3077  for (unsigned int d = 0; d < spacedim; ++d)
3078  vertices[i + v_offset][d] = mVertices[i][d];
3079 
3080  unsigned int valid_cell = c_offset;
3081  for (unsigned int i = 0; i < n_faces; ++i)
3082  {
3083  if (mFaces[i].mNumIndices == GeometryInfo<dim>::vertices_per_cell)
3084  {
3085  for (const unsigned int f : GeometryInfo<dim>::vertex_indices())
3086  {
3087  cells[valid_cell].vertices[f] =
3088  mFaces[i].mIndices[f] + v_offset;
3089  }
3090  cells[valid_cell].material_id = m;
3091  ++valid_cell;
3092  }
3093  else
3094  {
3095  AssertThrow(ignore_unsupported_types,
3096  ExcMessage("Face " + std::to_string(i) + " of mesh " +
3097  std::to_string(m) + " has " +
3098  std::to_string(mFaces[i].mNumIndices) +
3099  " vertices. We expected only " +
3102  }
3103  }
3104  cells.resize(valid_cell);
3105 
3106  // The vertices are added all at once. Cells are checked for
3107  // validity, so only valid_cells are now present in the deal.II
3108  // list of cells.
3109  v_offset += n_vertices;
3110  c_offset = valid_cell;
3111  }
3112 
3113  // No cells were read
3114  if (cells.size() == 0)
3115  return;
3116 
3117  if (remove_duplicates)
3118  {
3119  // The function delete_duplicated_vertices() needs to be called more
3120  // than once if a vertex is duplicated more than once. So we keep
3121  // calling it until the number of vertices does not change any more.
3122  unsigned int n_verts = 0;
3123  while (n_verts != vertices.size())
3124  {
3125  n_verts = vertices.size();
3126  std::vector<unsigned int> considered_vertices;
3128  vertices, cells, subcelldata, considered_vertices, tol);
3129  }
3130  }
3131 
3132  GridTools::delete_unused_vertices(vertices, cells, subcelldata);
3133  if (dim == spacedim)
3135  cells);
3136 
3138  if (dim == 2)
3139  tria->create_triangulation_compatibility(vertices, cells, subcelldata);
3140  else
3141  tria->create_triangulation(vertices, cells, subcelldata);
3142 #else
3143  (void)filename;
3144  (void)mesh_index;
3145  (void)remove_duplicates;
3146  (void)tol;
3147  (void)ignore_unsupported_types;
3148  AssertThrow(false, ExcNeedsAssimp());
3149 #endif
3150 }
3151 
3152 
3153 template <int dim, int spacedim>
3154 void
3156 {
3157  std::string line;
3158  while (in)
3159  {
3160  // get line
3161  getline(in, line);
3162 
3163  // check if this is a line that
3164  // consists only of spaces, and
3165  // if not put the whole thing
3166  // back and return
3167  if (std::find_if(line.begin(), line.end(), [](const char c) {
3168  return c != ' ';
3169  }) != line.end())
3170  {
3171  in.putback('\n');
3172  for (int i = line.length() - 1; i >= 0; --i)
3173  in.putback(line[i]);
3174  return;
3175  }
3176 
3177  // else: go on with next line
3178  }
3179 }
3180 
3181 
3182 
3183 template <int dim, int spacedim>
3184 void
3186  const char comment_start)
3187 {
3188  char c;
3189  // loop over the following comment
3190  // lines
3191  while (in.get(c) && c == comment_start)
3192  // loop over the characters after
3193  // the comment starter
3194  while (in.get() != '\n')
3195  ;
3196 
3197 
3198  // put back first character of
3199  // first non-comment line
3200  if (in)
3201  in.putback(c);
3202 
3203  // at last: skip additional empty lines, if present
3204  skip_empty_lines(in);
3205 }
3206 
3207 
3208 
3209 template <int dim, int spacedim>
3210 void
3212  const std::vector<CellData<dim>> & /*cells*/,
3213  const std::vector<Point<spacedim>> & /*vertices*/,
3214  std::ostream & /*out*/)
3215 {
3216  Assert(false, ExcNotImplemented());
3217 }
3218 
3219 
3220 
3221 template <>
3222 void
3223 GridIn<2>::debug_output_grid(const std::vector<CellData<2>> &cells,
3224  const std::vector<Point<2>> & vertices,
3225  std::ostream & out)
3226 {
3227  double min_x = vertices[cells[0].vertices[0]](0),
3228  max_x = vertices[cells[0].vertices[0]](0),
3229  min_y = vertices[cells[0].vertices[0]](1),
3230  max_y = vertices[cells[0].vertices[0]](1);
3231 
3232  for (unsigned int i = 0; i < cells.size(); ++i)
3233  {
3234  for (const auto vertex : cells[i].vertices)
3235  {
3236  const Point<2> &p = vertices[vertex];
3237 
3238  if (p(0) < min_x)
3239  min_x = p(0);
3240  if (p(0) > max_x)
3241  max_x = p(0);
3242  if (p(1) < min_y)
3243  min_y = p(1);
3244  if (p(1) > max_y)
3245  max_y = p(1);
3246  }
3247 
3248  out << "# cell " << i << std::endl;
3249  Point<2> center;
3250  for (const auto vertex : cells[i].vertices)
3251  center += vertices[vertex];
3252  center /= 4;
3253 
3254  out << "set label \"" << i << "\" at " << center(0) << ',' << center(1)
3255  << " center" << std::endl;
3256 
3257  // first two line right direction
3258  for (unsigned int f = 0; f < 2; ++f)
3259  out << "set arrow from " << vertices[cells[i].vertices[f]](0) << ','
3260  << vertices[cells[i].vertices[f]](1) << " to "
3261  << vertices[cells[i].vertices[(f + 1) % 4]](0) << ','
3262  << vertices[cells[i].vertices[(f + 1) % 4]](1) << std::endl;
3263  // other two lines reverse direction
3264  for (unsigned int f = 2; f < 4; ++f)
3265  out << "set arrow from " << vertices[cells[i].vertices[(f + 1) % 4]](0)
3266  << ',' << vertices[cells[i].vertices[(f + 1) % 4]](1) << " to "
3267  << vertices[cells[i].vertices[f]](0) << ','
3268  << vertices[cells[i].vertices[f]](1) << std::endl;
3269  out << std::endl;
3270  }
3271 
3272 
3273  out << std::endl
3274  << "set nokey" << std::endl
3275  << "pl [" << min_x << ':' << max_x << "][" << min_y << ':' << max_y
3276  << "] " << min_y << std::endl
3277  << "pause -1" << std::endl;
3278 }
3279 
3280 
3281 
3282 template <>
3283 void
3284 GridIn<3>::debug_output_grid(const std::vector<CellData<3>> &cells,
3285  const std::vector<Point<3>> & vertices,
3286  std::ostream & out)
3287 {
3288  for (const auto &cell : cells)
3289  {
3290  // line 0
3291  out << vertices[cell.vertices[0]] << std::endl
3292  << vertices[cell.vertices[1]] << std::endl
3293  << std::endl
3294  << std::endl;
3295  // line 1
3296  out << vertices[cell.vertices[1]] << std::endl
3297  << vertices[cell.vertices[2]] << std::endl
3298  << std::endl
3299  << std::endl;
3300  // line 2
3301  out << vertices[cell.vertices[3]] << std::endl
3302  << vertices[cell.vertices[2]] << std::endl
3303  << std::endl
3304  << std::endl;
3305  // line 3
3306  out << vertices[cell.vertices[0]] << std::endl
3307  << vertices[cell.vertices[3]] << std::endl
3308  << std::endl
3309  << std::endl;
3310  // line 4
3311  out << vertices[cell.vertices[4]] << std::endl
3312  << vertices[cell.vertices[5]] << std::endl
3313  << std::endl
3314  << std::endl;
3315  // line 5
3316  out << vertices[cell.vertices[5]] << std::endl
3317  << vertices[cell.vertices[6]] << std::endl
3318  << std::endl
3319  << std::endl;
3320  // line 6
3321  out << vertices[cell.vertices[7]] << std::endl
3322  << vertices[cell.vertices[6]] << std::endl
3323  << std::endl
3324  << std::endl;
3325  // line 7
3326  out << vertices[cell.vertices[4]] << std::endl
3327  << vertices[cell.vertices[7]] << std::endl
3328  << std::endl
3329  << std::endl;
3330  // line 8
3331  out << vertices[cell.vertices[0]] << std::endl
3332  << vertices[cell.vertices[4]] << std::endl
3333  << std::endl
3334  << std::endl;
3335  // line 9
3336  out << vertices[cell.vertices[1]] << std::endl
3337  << vertices[cell.vertices[5]] << std::endl
3338  << std::endl
3339  << std::endl;
3340  // line 10
3341  out << vertices[cell.vertices[2]] << std::endl
3342  << vertices[cell.vertices[6]] << std::endl
3343  << std::endl
3344  << std::endl;
3345  // line 11
3346  out << vertices[cell.vertices[3]] << std::endl
3347  << vertices[cell.vertices[7]] << std::endl
3348  << std::endl
3349  << std::endl;
3350  }
3351 }
3352 
3353 
3354 
3355 template <int dim, int spacedim>
3356 void
3357 GridIn<dim, spacedim>::read(const std::string &filename, Format format)
3358 {
3359  // Search file class for meshes
3360  PathSearch search("MESH");
3361  std::string name;
3362  // Open the file and remember its name
3363  if (format == Default)
3364  name = search.find(filename);
3365  else
3366  name = search.find(filename, default_suffix(format));
3367 
3368  std::ifstream in(name.c_str());
3369 
3370  if (format == Default)
3371  {
3372  const std::string::size_type slashpos = name.find_last_of('/');
3373  const std::string::size_type dotpos = name.find_last_of('.');
3374  if (dotpos < name.length() &&
3375  (dotpos > slashpos || slashpos == std::string::npos))
3376  {
3377  std::string ext = name.substr(dotpos + 1);
3378  format = parse_format(ext);
3379  }
3380  }
3381  if (format == netcdf)
3382  read_netcdf(filename);
3383  else
3384  read(in, format);
3385 }
3386 
3387 
3388 template <int dim, int spacedim>
3389 void
3390 GridIn<dim, spacedim>::read(std::istream &in, Format format)
3391 {
3392  if (format == Default)
3393  format = default_format;
3394 
3395  switch (format)
3396  {
3397  case dbmesh:
3398  read_dbmesh(in);
3399  return;
3400 
3401  case msh:
3402  read_msh(in);
3403  return;
3404 
3405  case vtk:
3406  read_vtk(in);
3407  return;
3408 
3409  case vtu:
3410  read_vtu(in);
3411  return;
3412 
3413  case unv:
3414  read_unv(in);
3415  return;
3416 
3417  case ucd:
3418  read_ucd(in);
3419  return;
3420 
3421  case abaqus:
3422  read_abaqus(in);
3423  return;
3424 
3425  case xda:
3426  read_xda(in);
3427  return;
3428 
3429  case netcdf:
3430  Assert(false,
3431  ExcMessage("There is no read_netcdf(istream &) function. "
3432  "Use the read_netcdf(string &filename) "
3433  "functions, instead."));
3434  return;
3435 
3436  case tecplot:
3437  read_tecplot(in);
3438  return;
3439 
3440  case assimp:
3441  Assert(false,
3442  ExcMessage("There is no read_assimp(istream &) function. "
3443  "Use the read_assimp(string &filename, ...) "
3444  "functions, instead."));
3445  return;
3446 
3447  case Default:
3448  break;
3449  }
3450  Assert(false, ExcInternalError());
3451 }
3452 
3453 
3454 
3455 template <int dim, int spacedim>
3456 std::string
3458 {
3459  switch (format)
3460  {
3461  case dbmesh:
3462  return ".dbmesh";
3463  case msh:
3464  return ".msh";
3465  case vtk:
3466  return ".vtk";
3467  case vtu:
3468  return ".vtu";
3469  case unv:
3470  return ".unv";
3471  case ucd:
3472  return ".inp";
3473  case abaqus:
3474  return ".inp"; // Typical suffix for Abaqus mesh files conflicts with
3475  // UCD.
3476  case xda:
3477  return ".xda";
3478  case netcdf:
3479  return ".nc";
3480  case tecplot:
3481  return ".dat";
3482  default:
3483  Assert(false, ExcNotImplemented());
3484  return ".unknown_format";
3485  }
3486 }
3487 
3488 
3489 
3490 template <int dim, int spacedim>
3492 GridIn<dim, spacedim>::parse_format(const std::string &format_name)
3493 {
3494  if (format_name == "dbmesh")
3495  return dbmesh;
3496 
3497  if (format_name == "msh")
3498  return msh;
3499 
3500  if (format_name == "unv")
3501  return unv;
3502 
3503  if (format_name == "vtk")
3504  return vtk;
3505 
3506  if (format_name == "vtu")
3507  return vtu;
3508 
3509  // This is also the typical extension of Abaqus input files.
3510  if (format_name == "inp")
3511  return ucd;
3512 
3513  if (format_name == "ucd")
3514  return ucd;
3515 
3516  if (format_name == "xda")
3517  return xda;
3518 
3519  if (format_name == "netcdf")
3520  return netcdf;
3521 
3522  if (format_name == "nc")
3523  return netcdf;
3524 
3525  if (format_name == "tecplot")
3526  return tecplot;
3527 
3528  if (format_name == "dat")
3529  return tecplot;
3530 
3531  if (format_name == "plt")
3532  // Actually, this is the extension for the
3533  // tecplot binary format, which we do not
3534  // support right now. However, some people
3535  // tend to create tecplot ascii files with
3536  // the extension 'plt' instead of
3537  // 'dat'. Thus, include this extension
3538  // here. If it actually is a binary file,
3539  // the read_tecplot() function will fail
3540  // and throw an exception, anyway.
3541  return tecplot;
3542 
3543  AssertThrow(false, ExcInvalidState());
3544  // return something weird
3545  return Format(Default);
3546 }
3547 
3548 
3549 
3550 template <int dim, int spacedim>
3551 std::string
3553 {
3554  return "dbmesh|msh|unv|vtk|vtu|ucd|abaqus|xda|netcdf|tecplot|assimp";
3555 }
3556 
3557 
3558 
3559 namespace
3560 {
3561  template <int dim, int spacedim>
3562  Abaqus_to_UCD<dim, spacedim>::Abaqus_to_UCD()
3563  : tolerance(5e-16) // Used to offset Cubit tolerance error when outputting
3564  // value close to zero
3565  {
3566  AssertThrow(spacedim == 2 || spacedim == 3, ExcNotImplemented());
3567  }
3568 
3569 
3570 
3571  // Convert from a string to some other data type
3572  // Reference: http://www.codeguru.com/forum/showthread.php?t=231054
3573  template <class T>
3574  bool
3575  from_string(T &t, const std::string &s, std::ios_base &(*f)(std::ios_base &))
3576  {
3577  std::istringstream iss(s);
3578  return !(iss >> f >> t).fail();
3579  }
3580 
3581 
3582 
3583  // Extract an integer from a string
3584  int
3585  extract_int(const std::string &s)
3586  {
3587  std::string tmp;
3588  for (const char c : s)
3589  {
3590  if (isdigit(c))
3591  {
3592  tmp += c;
3593  }
3594  }
3595 
3596  int number = 0;
3597  from_string(number, tmp, std::dec);
3598  return number;
3599  }
3600 
3601 
3602 
3603  template <int dim, int spacedim>
3604  void
3605  Abaqus_to_UCD<dim, spacedim>::read_in_abaqus(std::istream &input_stream)
3606  {
3607  // References:
3608  // http://www.egr.msu.edu/software/abaqus/Documentation/docs/v6.7/books/usb/default.htm?startat=pt01ch02.html
3609  // http://www.cprogramming.com/tutorial/string.html
3610 
3611  AssertThrow(input_stream, ExcIO());
3612  std::string line;
3613 
3614  while (std::getline(input_stream, line))
3615  {
3616  cont:
3617  std::transform(line.begin(), line.end(), line.begin(), ::toupper);
3618 
3619  if (line.compare("*HEADING") == 0 || line.compare(0, 2, "**") == 0 ||
3620  line.compare(0, 5, "*PART") == 0)
3621  {
3622  // Skip header and comments
3623  while (std::getline(input_stream, line))
3624  {
3625  if (line[0] == '*')
3626  goto cont; // My eyes, they burn!
3627  }
3628  }
3629  else if (line.compare(0, 5, "*NODE") == 0)
3630  {
3631  // Extract list of vertices
3632  // Header line might be:
3633  // *NODE, NSET=ALLNODES
3634  // *NODE
3635 
3636  // Contains lines in the form:
3637  // Index, x, y, z
3638  while (std::getline(input_stream, line))
3639  {
3640  if (line[0] == '*')
3641  goto cont;
3642 
3643  std::vector<double> node(spacedim + 1);
3644 
3645  std::istringstream iss(line);
3646  char comma;
3647  for (unsigned int i = 0; i < spacedim + 1; ++i)
3648  iss >> node[i] >> comma;
3649 
3650  node_list.push_back(node);
3651  }
3652  }
3653  else if (line.compare(0, 8, "*ELEMENT") == 0)
3654  {
3655  // Element construction.
3656  // There are different header formats, the details
3657  // of which we're not particularly interested in except
3658  // whether they represent quads or hexahedrals.
3659  // *ELEMENT, TYPE=S4R, ELSET=EB<material id>
3660  // *ELEMENT, TYPE=C3D8R, ELSET=EB<material id>
3661  // *ELEMENT, TYPE=C3D8
3662  // Elements itself (n=4 or n=8):
3663  // Index, i[0], ..., i[n]
3664 
3665  int material = 0;
3666  // Scan for material id
3667  {
3668  const std::string before_material = "ELSET=EB";
3669  const std::size_t idx = line.find(before_material);
3670  if (idx != std::string::npos)
3671  {
3672  from_string(material,
3673  line.substr(idx + before_material.size()),
3674  std::dec);
3675  }
3676  }
3677 
3678  // Read ELEMENT definition
3679  while (std::getline(input_stream, line))
3680  {
3681  if (line[0] == '*')
3682  goto cont;
3683 
3684  std::istringstream iss(line);
3685  char comma;
3686 
3687  // We will store the material id in the zeroth entry of the
3688  // vector and the rest of the elements represent the global
3689  // node numbers
3690  const unsigned int n_data_per_cell =
3692  std::vector<double> cell(n_data_per_cell);
3693  for (unsigned int i = 0; i < n_data_per_cell; ++i)
3694  iss >> cell[i] >> comma;
3695 
3696  // Overwrite cell index from file by material
3697  cell[0] = static_cast<double>(material);
3698  cell_list.push_back(cell);
3699  }
3700  }
3701  else if (line.compare(0, 8, "*SURFACE") == 0)
3702  {
3703  // Extract the definitions of boundary surfaces
3704  // Old format from Cubit:
3705  // *SURFACE, NAME=SS<boundary indicator>
3706  // <element index>, S<face number>
3707  // Abaqus default format:
3708  // *SURFACE, TYPE=ELEMENT, NAME=SURF-<indicator>
3709 
3710  // Get name of the surface and extract id from it;
3711  // this will be the boundary indicator
3712  const std::string name_key = "NAME=";
3713  const std::size_t name_idx_start =
3714  line.find(name_key) + name_key.size();
3715  std::size_t name_idx_end = line.find(',', name_idx_start);
3716  if (name_idx_end == std::string::npos)
3717  {
3718  name_idx_end = line.size();
3719  }
3720  const int b_indicator = extract_int(
3721  line.substr(name_idx_start, name_idx_end - name_idx_start));
3722 
3723  // Read SURFACE definition
3724  // Note that the orientation of the faces is embedded within the
3725  // definition of each "set" of faces that comprise the surface
3726  // These are either marked by an "S" or "E" in 3d or 2d
3727  // respectively.
3728  while (std::getline(input_stream, line))
3729  {
3730  if (line[0] == '*')
3731  goto cont;
3732 
3733  // Change all characters to upper case
3734  std::transform(line.begin(),
3735  line.end(),
3736  line.begin(),
3737  ::toupper);
3738 
3739  // Surface can be created from ELSET, or directly from cells
3740  // If elsets_list contains a key with specific name - refers to
3741  // that ELSET, otherwise refers to cell
3742  std::istringstream iss(line);
3743  int el_idx;
3744  int face_number;
3745  char temp;
3746 
3747  // Get relevant faces, taking into account the element
3748  // orientation
3749  std::vector<double> quad_node_list;
3750  const std::string elset_name = line.substr(0, line.find(','));
3751  if (elsets_list.count(elset_name) != 0)
3752  {
3753  // Surface refers to ELSET
3754  std::string stmp;
3755  iss >> stmp >> temp >> face_number;
3756 
3757  const std::vector<int> cells = elsets_list[elset_name];
3758  for (const int cell : cells)
3759  {
3760  el_idx = cell;
3761  quad_node_list =
3762  get_global_node_numbers(el_idx, face_number);
3763  quad_node_list.insert(quad_node_list.begin(),
3764  b_indicator);
3765 
3766  face_list.push_back(quad_node_list);
3767  }
3768  }
3769  else
3770  {
3771  // Surface refers directly to elements
3772  char comma;
3773  iss >> el_idx >> comma >> temp >> face_number;
3774  quad_node_list =
3775  get_global_node_numbers(el_idx, face_number);
3776  quad_node_list.insert(quad_node_list.begin(), b_indicator);
3777 
3778  face_list.push_back(quad_node_list);
3779  }
3780  }
3781  }
3782  else if (line.compare(0, 6, "*ELSET") == 0)
3783  {
3784  // Get ELSET name.
3785  // Materials are attached to elsets with specific name
3786  std::string elset_name;
3787  {
3788  const std::string elset_key = "*ELSET, ELSET=";
3789  const std::size_t idx = line.find(elset_key);
3790  if (idx != std::string::npos)
3791  {
3792  const std::string comma = ",";
3793  const std::size_t first_comma = line.find(comma);
3794  const std::size_t second_comma =
3795  line.find(comma, first_comma + 1);
3796  const std::size_t elset_name_start =
3797  line.find(elset_key) + elset_key.size();
3798  elset_name = line.substr(elset_name_start,
3799  second_comma - elset_name_start);
3800  }
3801  }
3802 
3803  // There are two possibilities of storing cells numbers in ELSET:
3804  // 1. If the header contains the 'GENERATE' keyword, then the next
3805  // line describes range of cells as:
3806  // cell_id_start, cell_id_end, cell_step
3807  // 2. If the header does not contain the 'GENERATE' keyword, then
3808  // the next lines contain cells numbers
3809  std::vector<int> elements;
3810  const std::size_t generate_idx = line.find("GENERATE");
3811  if (generate_idx != std::string::npos)
3812  {
3813  // Option (1)
3814  std::getline(input_stream, line);
3815  std::istringstream iss(line);
3816  char comma;
3817  int elid_start;
3818  int elid_end;
3819  int elis_step = 1; // Default if case stride not provided
3820 
3821  // Some files don't have the stride size
3822  // Compare mesh test cases ./grids/abaqus/3d/other_simple.inp to
3823  // ./grids/abaqus/2d/2d_test_abaqus.inp
3824  iss >> elid_start >> comma >> elid_end;
3825  AssertThrow(comma == ',',
3826  ExcMessage(
3827  std::string(
3828  "While reading an ABAQUS file, the reader "
3829  "expected a comma but found a <") +
3830  comma + "> in the line <" + line + ">."));
3831  AssertThrow(
3832  elid_start <= elid_end,
3833  ExcMessage(
3834  std::string(
3835  "While reading an ABAQUS file, the reader encountered "
3836  "a GENERATE statement in which the upper bound <") +
3837  Utilities::int_to_string(elid_end) +
3838  "> for the element numbers is not larger or equal "
3839  "than the lower bound <" +
3840  Utilities::int_to_string(elid_start) + ">."));
3841 
3842  // https://stackoverflow.com/questions/8046357/how-do-i-check-if-a-stringstream-variable-is-empty-null
3843  if (iss.rdbuf()->in_avail() != 0)
3844  iss >> comma >> elis_step;
3845  AssertThrow(comma == ',',
3846  ExcMessage(
3847  std::string(
3848  "While reading an ABAQUS file, the reader "
3849  "expected a comma but found a <") +
3850  comma + "> in the line <" + line + ">."));
3851 
3852  for (int i = elid_start; i <= elid_end; i += elis_step)
3853  elements.push_back(i);
3854  elsets_list[elset_name] = elements;
3855 
3856  std::getline(input_stream, line);
3857  }
3858  else
3859  {
3860  // Option (2)
3861  while (std::getline(input_stream, line))
3862  {
3863  if (line[0] == '*')
3864  break;
3865 
3866  std::istringstream iss(line);
3867  char comma;
3868  int elid;
3869  while (!iss.eof())
3870  {
3871  iss >> elid >> comma;
3872  AssertThrow(
3873  comma == ',',
3874  ExcMessage(
3875  std::string(
3876  "While reading an ABAQUS file, the reader "
3877  "expected a comma but found a <") +
3878  comma + "> in the line <" + line + ">."));
3879 
3880  elements.push_back(elid);
3881  }
3882  }
3883 
3884  elsets_list[elset_name] = elements;
3885  }
3886 
3887  goto cont;
3888  }
3889  else if (line.compare(0, 5, "*NSET") == 0)
3890  {
3891  // Skip nodesets; we have no use for them
3892  while (std::getline(input_stream, line))
3893  {
3894  if (line[0] == '*')
3895  goto cont;
3896  }
3897  }
3898  else if (line.compare(0, 14, "*SOLID SECTION") == 0)
3899  {
3900  // The ELSET name, which describes a section for particular material
3901  const std::string elset_key = "ELSET=";
3902  const std::size_t elset_start =
3903  line.find("ELSET=") + elset_key.size();
3904  const std::size_t elset_end = line.find(',', elset_start + 1);
3905  const std::string elset_name =
3906  line.substr(elset_start, elset_end - elset_start);
3907 
3908  // Solid material definition.
3909  // We assume that material id is taken from material name,
3910  // eg. "Material-1" -> ID=1
3911  const std::string material_key = "MATERIAL=";
3912  const std::size_t last_equal =
3913  line.find("MATERIAL=") + material_key.size();
3914  const std::size_t material_id_start = line.find('-', last_equal);
3915  int material_id = 0;
3916  from_string(material_id,
3917  line.substr(material_id_start + 1),
3918  std::dec);
3919 
3920  // Assign material id to cells
3921  const std::vector<int> &elset_cells = elsets_list[elset_name];
3922  for (const int elset_cell : elset_cells)
3923  {
3924  const int cell_id = elset_cell - 1;
3925  cell_list[cell_id][0] = material_id;
3926  }
3927  }
3928  // Note: All other lines / entries are ignored
3929  }
3930  }
3931 
3932  template <int dim, int spacedim>
3933  std::vector<double>
3934  Abaqus_to_UCD<dim, spacedim>::get_global_node_numbers(
3935  const int face_cell_no,
3936  const int face_cell_face_no) const
3937  {
3938  std::vector<double> quad_node_list(GeometryInfo<dim>::vertices_per_face);
3939 
3940  // These orderings were reverse engineered by hand and may
3941  // conceivably be erroneous.
3942  // TODO: Currently one test (2d unstructured mesh) in the test
3943  // suite fails, presumably because of an ordering issue.
3944  if (dim == 2)
3945  {
3946  if (face_cell_face_no == 1)
3947  {
3948  quad_node_list[0] = cell_list[face_cell_no - 1][1];
3949  quad_node_list[1] = cell_list[face_cell_no - 1][2];
3950  }
3951  else if (face_cell_face_no == 2)
3952  {
3953  quad_node_list[0] = cell_list[face_cell_no - 1][2];
3954  quad_node_list[1] = cell_list[face_cell_no - 1][3];
3955  }
3956  else if (face_cell_face_no == 3)
3957  {
3958  quad_node_list[0] = cell_list[face_cell_no - 1][3];
3959  quad_node_list[1] = cell_list[face_cell_no - 1][4];
3960  }
3961  else if (face_cell_face_no == 4)
3962  {
3963  quad_node_list[0] = cell_list[face_cell_no - 1][4];
3964  quad_node_list[1] = cell_list[face_cell_no - 1][1];
3965  }
3966  else
3967  {
3968  AssertThrow(face_cell_face_no <= 4,
3969  ExcMessage("Invalid face number in 2d"));
3970  }
3971  }
3972  else if (dim == 3)
3973  {
3974  if (face_cell_face_no == 1)
3975  {
3976  quad_node_list[0] = cell_list[face_cell_no - 1][1];
3977  quad_node_list[1] = cell_list[face_cell_no - 1][4];
3978  quad_node_list[2] = cell_list[face_cell_no - 1][3];
3979  quad_node_list[3] = cell_list[face_cell_no - 1][2];
3980  }
3981  else if (face_cell_face_no == 2)
3982  {
3983  quad_node_list[0] = cell_list[face_cell_no - 1][5];
3984  quad_node_list[1] = cell_list[face_cell_no - 1][8];
3985  quad_node_list[2] = cell_list[face_cell_no - 1][7];
3986  quad_node_list[3] = cell_list[face_cell_no - 1][6];
3987  }
3988  else if (face_cell_face_no == 3)
3989  {
3990  quad_node_list[0] = cell_list[face_cell_no - 1][1];
3991  quad_node_list[1] = cell_list[face_cell_no - 1][2];
3992  quad_node_list[2] = cell_list[face_cell_no - 1][6];
3993  quad_node_list[3] = cell_list[face_cell_no - 1][5];
3994  }
3995  else if (face_cell_face_no == 4)
3996  {
3997  quad_node_list[0] = cell_list[face_cell_no - 1][2];
3998  quad_node_list[1] = cell_list[face_cell_no - 1][3];
3999  quad_node_list[2] = cell_list[face_cell_no - 1][7];
4000  quad_node_list[3] = cell_list[face_cell_no - 1][6];
4001  }
4002  else if (face_cell_face_no == 5)
4003  {
4004  quad_node_list[0] = cell_list[face_cell_no - 1][3];
4005  quad_node_list[1] = cell_list[face_cell_no - 1][4];
4006  quad_node_list[2] = cell_list[face_cell_no - 1][8];
4007  quad_node_list[3] = cell_list[face_cell_no - 1][7];
4008  }
4009  else if (face_cell_face_no == 6)
4010  {
4011  quad_node_list[0] = cell_list[face_cell_no - 1][1];
4012  quad_node_list[1] = cell_list[face_cell_no - 1][5];
4013  quad_node_list[2] = cell_list[face_cell_no - 1][8];
4014  quad_node_list[3] = cell_list[face_cell_no - 1][4];
4015  }
4016  else
4017  {
4018  AssertThrow(face_cell_no <= 6,
4019  ExcMessage("Invalid face number in 3d"));
4020  }
4021  }
4022  else
4023  {
4024  AssertThrow(dim == 2 || dim == 3, ExcNotImplemented());
4025  }
4026 
4027  return quad_node_list;
4028  }
4029 
4030  template <int dim, int spacedim>
4031  void
4032  Abaqus_to_UCD<dim, spacedim>::write_out_avs_ucd(std::ostream &output) const
4033  {
4034  // References:
4035  // http://www.dealii.org/developer/doxygen/deal.II/structGeometryInfo.html
4036  // http://people.scs.fsu.edu/~burkardt/data/ucd/ucd.html
4037 
4038  AssertThrow(output, ExcIO());
4039 
4040  // save old formatting options
4041  const boost::io::ios_base_all_saver formatting_saver(output);
4042 
4043  // Write out title - Note: No other commented text can be inserted below the
4044  // title in a UCD file
4045  output << "# Abaqus to UCD mesh conversion" << std::endl;
4046  output << "# Mesh type: AVS UCD" << std::endl;
4047 
4048  // ========================================================
4049  // ASCII UCD File Format
4050  // The input file cannot contain blank lines or lines with leading blanks.
4051  // Comments, if present, must precede all data in the file.
4052  // Comments within the data will cause read errors.
4053  // The general order of the data is as follows:
4054  // 1. Numbers defining the overall structure, including the number of nodes,
4055  // the number of cells, and the length of the vector of data associated
4056  // with the nodes, cells, and the model.
4057  // e.g. 1:
4058  // <num_nodes> <num_cells> <num_ndata> <num_cdata> <num_mdata>
4059  // e.g. 2:
4060  // n_elements = n_hex_cells + n_bc_quads + n_quad_cells + n_bc_edges
4061  // outfile.write(str(n_nodes) + " " + str(n_elements) + " 0 0 0\n")
4062  // 2. For each node, its node id and the coordinates of that node in space.
4063  // Node-ids must be integers, but any number including non sequential
4064  // numbers can be used. Mid-edge nodes are treated like any other node.
4065  // 3. For each cell: its cell-id, material, cell type (hexahedral, pyramid,
4066  // etc.), and the list of node-ids that correspond to each of the cell's
4067  // vertices. The below table specifies the different cell types and the
4068  // keyword used to represent them in the file.
4069 
4070  // Write out header
4071  output << node_list.size() << "\t" << (cell_list.size() + face_list.size())
4072  << "\t0\t0\t0" << std::endl;
4073 
4074  output.width(16);
4075  output.precision(8);
4076 
4077  // Write out node numbers
4078  // Loop over all nodes
4079  for (const auto &node : node_list)
4080  {
4081  // Node number
4082  output << node[0] << "\t";
4083 
4084  // Node coordinates
4085  output.setf(std::ios::scientific, std::ios::floatfield);
4086  for (unsigned int jj = 1; jj < spacedim + 1; ++jj)
4087  {
4088  // invoke tolerance -> set points close to zero equal to zero
4089  if (std::abs(node[jj]) > tolerance)
4090  output << static_cast<double>(node[jj]) << "\t";
4091  else
4092  output << 0.0 << "\t";
4093  }
4094  if (spacedim == 2)
4095  output << 0.0 << "\t";
4096 
4097  output << std::endl;
4098  output.unsetf(std::ios::floatfield);
4099  }
4100 
4101  // Write out cell node numbers
4102  for (unsigned int ii = 0; ii < cell_list.size(); ++ii)
4103  {
4104  output << ii + 1 << "\t" << cell_list[ii][0] << "\t"
4105  << (dim == 2 ? "quad" : "hex") << "\t";
4106  for (unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_cell + 1;
4107  ++jj)
4108  output << cell_list[ii][jj] << "\t";
4109 
4110  output << std::endl;
4111  }
4112 
4113  // Write out quad node numbers
4114  for (unsigned int ii = 0; ii < face_list.size(); ++ii)
4115  {
4116  output << ii + 1 << "\t" << face_list[ii][0] << "\t"
4117  << (dim == 2 ? "line" : "quad") << "\t";
4118  for (unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_face + 1;
4119  ++jj)
4120  output << face_list[ii][jj] << "\t";
4121 
4122  output << std::endl;
4123  }
4124  }
4125 } // namespace
4126 
4127 
4128 // explicit instantiations
4129 #include "grid_in.inst"
4130 
GridIn::default_suffix
static std::string default_suffix(const Format format)
Definition: grid_in.cc:3457
GridIn::debug_output_grid
static void debug_output_grid(const std::vector< CellData< dim >> &cells, const std::vector< Point< spacedim >> &vertices, std::ostream &out)
Definition: grid_in.cc:3211
GridIn::read
void read(std::istream &in, Format format=Default)
Definition: grid_in.cc:3390
StandardExceptions::ExcNeedsAssimp
static ::ExceptionBase & ExcNeedsAssimp()
GridIn::read_ucd
void read_ucd(std::istream &in, const bool apply_all_indicators_to_manifolds=false)
Definition: grid_in.cc:781
StandardExceptions::ExcIO
static ::ExceptionBase & ExcIO()
GridIn::read_xda
void read_xda(std::istream &in)
Definition: grid_in.cc:1269
SubCellData::boundary_quads
std::vector< CellData< 2 > > boundary_quads
Definition: tria_description.h:215
GridIn::read_abaqus
void read_abaqus(std::istream &in, const bool apply_all_indicators_to_manifolds=false)
Definition: grid_in.cc:1046
Utilities::int_to_string
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:474
PathSearch
Definition: path_search.h:84
StandardExceptions::ExcNeedsNetCDF
static ::ExceptionBase & ExcNeedsNetCDF()
GridIn::attach_triangulation
void attach_triangulation(Triangulation< dim, spacedim > &tria)
Definition: grid_in.cc:106
CellData
Definition: tria_description.h:67
GridIn::parse_format
static Format parse_format(const std::string &format_name)
Definition: grid_in.cc:3492
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Triangulation< 1, spacedim >
GridTools::delete_unused_vertices
void delete_unused_vertices(std::vector< Point< spacedim >> &vertices, std::vector< CellData< dim >> &cells, SubCellData &subcelldata)
Definition: grid_tools.cc:513
tria.h
DEAL_II_FALLTHROUGH
#define DEAL_II_FALLTHROUGH
Definition: config.h:153
parallel::transform
void transform(const InputIterator &begin_in, const InputIterator &end_in, OutputIterator out, const Predicate &predicate, const unsigned int grainsize)
Definition: parallel.h:211
StandardExceptions::ExcInvalidState
static ::ExceptionBase & ExcInvalidState()
utilities.h
DataOutBase::vtk
@ vtk
Definition: data_out_base.h:1605
Utilities::decode_base64
std::vector< unsigned char > decode_base64(const std::string &base64_input)
Definition: utilities.cc:450
GeometryInfo
Definition: geometry_info.h:1224
DataOutBase::vtu
@ vtu
Definition: data_out_base.h:1610
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
DataOutBase::default_format
@ default_format
Definition: data_out_base.h:1552
Patterns::Tools::to_string
std::string to_string(const T &t)
Definition: patterns.h:2360
GridIn::read_dbmesh
void read_dbmesh(std::istream &in)
Definition: grid_in.cc:1102
Physics::Elasticity::Kinematics::e
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
numbers::flat_manifold_id
const types::manifold_id flat_manifold_id
Definition: types.h:273
GridReordering::reorder_cells
static void reorder_cells(std::vector< CellData< dim >> &original_cells, const bool use_new_style_ordering=false)
Definition: grid_reordering.cc:1078
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
internal::TriangulationImplementation::ExcInvalidVertexIndex
static ::ExceptionBase & ExcInvalidVertexIndex(int arg1, int arg2, int arg3)
grid_in.h
Utilities::get_integer_at_position
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
Definition: utilities.cc:854
GridIn::Format
Format
Definition: grid_in.h:314
GridIn::GridIn
GridIn()
Definition: grid_in.cc:98
GridIn::read_unv
void read_unv(std::istream &in)
Definition: grid_in.cc:532
GridTools::delete_duplicated_vertices
void delete_duplicated_vertices(std::vector< Point< spacedim >> &all_vertices, std::vector< CellData< dim >> &cells, SubCellData &subcelldata, std::vector< unsigned int > &considered_vertices, const double tol=1e-12)
Definition: grid_tools.cc:614
LAPACKSupport::T
static const char T
Definition: lapack_support.h:163
types::material_id
unsigned int material_id
Definition: types.h:152
StandardExceptions::ExcIndexRange
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
grid_tools.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
GridIn::read_vtu
void read_vtu(std::istream &in)
Definition: grid_in.cc:505
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
GridIn::get_format_names
static std::string get_format_names()
Definition: grid_in.cc:3552
numbers::internal_face_boundary_id
const types::boundary_id internal_face_boundary_id
Definition: types.h:250
GridIn::parse_tecplot_header
static void parse_tecplot_header(std::string &header, std::vector< unsigned int > &tecplot2deal, unsigned int &n_vars, unsigned int &n_vertices, unsigned int &n_cells, std::vector< unsigned int > &IJK, bool &structured, bool &blocked)
Definition: grid_in.cc:2514
exceptions.h
DataOutBase::tecplot
@ tecplot
Definition: data_out_base.h:1592
numbers::invalid_material_id
const types::material_id invalid_material_id
Definition: types.h:223
GridIn::read_assimp
void read_assimp(const std::string &filename, const unsigned int mesh_index=numbers::invalid_unsigned_int, const bool remove_duplicates=true, const double tol=1e-12, const bool ignore_unsupported_element_types=true)
Definition: grid_in.cc:2998
internal::dummy
const types::global_dof_index * dummy()
Definition: dof_handler.cc:59
unsigned int
vertices
Point< 3 > vertices[4]
Definition: data_out_base.cc:174
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
LinearAlgebra::CUDAWrappers::kernel::size_type
types::global_dof_index size_type
Definition: cuda_kernels.h:45
DataOutBase::default_suffix
std::string default_suffix(const OutputFormat output_format)
Definition: data_out_base.cc:2329
vertex_indices
std::vector< unsigned int > vertex_indices
Definition: tria.cc:2244
GridIn::skip_comment_lines
static void skip_comment_lines(std::istream &in, const char comment_start)
Definition: grid_in.cc:3185
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
grid_reordering.h
SubCellData::check_consistency
bool check_consistency(const unsigned int dim) const
Definition: tria.cc:82
StandardExceptions::ExcImpossibleInDim
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
GridIn::read_vtk
void read_vtk(std::istream &in)
Definition: grid_in.cc:115
path_search.h
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
GridIn::skip_empty_lines
static void skip_empty_lines(std::istream &in)
Definition: grid_in.cc:3155
Utilities::MPI::min
T min(const T &t, const MPI_Comm &mpi_communicator)
Point< 3 >
Exceptions::DataOutImplementation::ExcNoTriangulationSelected
static ::ExceptionBase & ExcNoTriangulationSelected()
GridIn::read_msh
void read_msh(std::istream &in)
Definition: grid_in.cc:1427
triangulation
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Definition: p4est_wrappers.cc:69
Utilities::match_at_string_start
bool match_at_string_start(const std::string &name, const std::string &pattern)
Definition: utilities.cc:839
SubCellData::boundary_lines
std::vector< CellData< 1 > > boundary_lines
Definition: tria_description.h:207
SubCellData
Definition: tria_description.h:199
internal::TriangulationImplementation::n_cells
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
Definition: tria.cc:12618
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Utilities::break_text_into_lines
std::vector< std::string > break_text_into_lines(const std::string &original_text, const unsigned int width, const char delimiter=' ')
Definition: utilities.cc:759
PathSearch::find
std::string find(const std::string &filename, const char *open_mode="r")
Definition: path_search.cc:171
DataOutBase::ucd
@ ucd
Definition: data_out_base.h:1567
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
center
Point< 3 > center
Definition: data_out_base.cc:169
GridIn::read_netcdf
void read_netcdf(const std::string &filename)
Utilities::decompress
std::string decompress(const std::string &compressed_input)
Definition: utilities.cc:415
GridIn::read_tecplot
void read_tecplot(std::istream &in)
Definition: grid_in.cc:2989
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
GridReordering::invert_all_cells_of_negative_grid
static void invert_all_cells_of_negative_grid(const std::vector< Point< spacedim >> &all_vertices, std::vector< CellData< dim >> &original_cells)