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