Reference documentation for deal.II version GIT relicensing-1321-g19b0133728 2024-07-26 07:50:01+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
grid_in.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1999 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_grid_in_h
16#define dealii_grid_in_h
17
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/point.h>
24
25#include <iostream>
26#include <string>
27#include <vector>
28
30
31// Forward declarations
32#ifndef DOXYGEN
33template <int dim, int spacedim>
35class Triangulation;
36template <int dim>
37struct CellData;
38#endif
39
309template <int dim, int spacedim = dim>
311{
312public:
344
348 GridIn();
349
354
358 void
360
365 void
366 read(std::istream &in, Format format = Default);
367
372 void
373 read(const std::string &in, Format format = Default);
374
415 void
416 read_vtk(std::istream &in);
417
432 void
433 read_vtu(std::istream &in);
434
435
443 void
444 read_unv(std::istream &in);
445
455 void
456 read_ucd(std::istream &in,
457 const bool apply_all_indicators_to_manifolds = false);
458
498 void
499 read_abaqus(std::istream &in,
500 const bool apply_all_indicators_to_manifolds = false);
501
505 void
506 read_dbmesh(std::istream &in);
507
511 void
512 read_xda(std::istream &in);
513
521 void
522 read_msh(std::istream &in);
523
524#ifdef DEAL_II_GMSH_WITH_API
582 void
583 read_msh(const std::string &filename);
584#endif
585
622 void
623 read_comsol_mphtxt(std::istream &in);
624
629 void
630 read_tecplot(std::istream &in);
631
666 void
667 read_assimp(const std::string &filename,
668 const unsigned int mesh_index = numbers::invalid_unsigned_int,
669 const bool remove_duplicates = true,
670 const double tol = 1e-12,
671 const bool ignore_unsupported_element_types = true);
672
682 {
687 std::vector<std::vector<int>> id_to_sideset_ids;
688 };
689
747 read_exodusii(const std::string &filename,
748 const bool apply_all_indicators_to_manifolds = false);
749
753 static std::string
754 default_suffix(const Format format);
755
759 static Format
760 parse_format(const std::string &format_name);
761
767 static std::string
769
774 int,
775 << "The section type <" << arg1 << "> in an UNV "
776 << "input file is not implemented.");
777
782 int,
783 << "The element type <" << arg1 << "> in an UNV "
784 << "input file is not implemented.");
785
790 std::string,
791 << "The identifier <" << arg1 << "> as name of a "
792 << "part in an UCD input file is unknown or the "
793 << "respective input routine is not implemented."
794 << "(Maybe the space dimension of triangulation and "
795 << "input file do not match?");
800 "No Triangulation has been attached to this GridIn object "
801 "so that nothing can be filled during any read function "
802 "calls. Please pass a reference to the Triangulation tria "
803 "to be filled in the constructor GridIn(tria) or attach "
804 "it with the function call GridIn::attach_triangulation().");
810 int,
811 int,
812 << "While creating cell " << arg1
813 << ", you are referencing a vertex with index " << arg2
814 << " but no vertex with this index has been described in the input file.");
820 int,
821 int,
822 int,
823 << "While creating cell " << arg1 << " (which is numbered as " << arg2
824 << " in the input file), you are referencing a vertex with index " << arg3
825 << " but no vertex with this index has been described in the input file.");
834 std::string,
835 << "The string <" << arg1
836 << "> is not recognized at the present"
837 << " position of a DB Mesh file.");
838
844 int,
845 << "The specified dimension " << arg1
846 << " is not the same as that of the triangulation to be created.");
847
849 std::string,
850 << "The string <" << arg1
851 << "> is not recognized at the present"
852 << " position of a Gmsh Mesh file.");
853
855 int,
856 << "The Element Identifier <" << arg1 << "> is not "
857 << "supported in the deal.II library when "
858 << "reading meshes in " << dim << " dimensions.\n"
859 << "Supported elements are: \n"
860 << "ELM-TYPE\n"
861 << "1 Line (2 nodes, 1 edge).\n"
862 << "2 Triangle (3 nodes, 3 edges).\n"
863 << "3 Quadrilateral (4 nodes, 4 edges).\n"
864 << "4 Tetrahedron (4 nodes, 6 edges, 4 faces) when in 3d.\n"
865 << "5 Hexahedron (8 nodes, 12 edges, 6 faces) when in 3d.\n"
866 << "15 Point (1 node, ignored when read).");
867
868
871 unsigned int,
872 unsigned int,
873 "While reading a gmsh file, the reader function did not find "
874 "any cells. This sometimes happens if the file only contains a "
875 "surface mesh, but not a volume mesh."
876 "\n\n"
877 "The reader function did find " +
878 std::to_string(arg1) + " lines and " + std::to_string(arg2) +
879 " facets (surface triangles or quadrilaterals).");
880
881protected:
886
907 static void
908 debug_output_grid(const std::vector<CellData<dim>> &cells,
909 const std::vector<Point<spacedim>> &vertices,
910 std::ostream &out);
911
912private:
917 static void
918 skip_empty_lines(std::istream &in);
919
927 static void
928 skip_comment_lines(std::istream &in, const char comment_start);
929
937 static void
938 parse_tecplot_header(std::string &header,
939 std::vector<unsigned int> &tecplot2deal,
940 unsigned int &n_vars,
941 unsigned int &n_vertices,
942 unsigned int &n_cells,
943 std::vector<unsigned int> &IJK,
944 bool &structured,
945 bool &blocked);
946
951};
952
953/* -------------- declaration of explicit specializations ------------- */
954
955#ifndef DOXYGEN
956
957template <>
958void
959GridIn<2>::debug_output_grid(const std::vector<CellData<2>> &cells,
960 const std::vector<Point<2>> &vertices,
961 std::ostream &out);
962
963
964template <>
965void
966GridIn<2, 3>::debug_output_grid(const std::vector<CellData<2>> &cells,
967 const std::vector<Point<3>> &vertices,
968 std::ostream &out);
969template <>
970void
971GridIn<3>::debug_output_grid(const std::vector<CellData<3>> &cells,
972 const std::vector<Point<3>> &vertices,
973 std::ostream &out);
974#endif // DOXYGEN
975
977
978#endif
@ Default
Use GridIn::default_format stored in this object.
Definition grid_in.h:320
@ vtu
Use read_vtu()
Definition grid_in.h:338
@ xda
Use read_xda()
Definition grid_in.h:330
@ tecplot
Use read_tecplot()
Definition grid_in.h:334
@ msh
Use read_msh()
Definition grid_in.h:332
@ ucd
Use read_ucd()
Definition grid_in.h:324
@ dbmesh
Use read_dbmesh()
Definition grid_in.h:328
@ abaqus
Use read_abaqus()
Definition grid_in.h:326
@ assimp
Use read_assimp()
Definition grid_in.h:340
@ unv
Use read_unv()
Definition grid_in.h:322
@ exodusii
Use read_exodusii()
Definition grid_in.h:342
@ vtk
Use read_vtk()
Definition grid_in.h:336
SmartPointer< Triangulation< dim, spacedim >, GridIn< dim, spacedim > > tria
Definition grid_in.h:885
void read_vtk(std::istream &in)
Definition grid_in.cc:162
GridIn()
Definition grid_in.cc:136
static void skip_empty_lines(std::istream &in)
Definition grid_in.cc:3937
void read_assimp(const std::string &filename, const unsigned int mesh_index=numbers::invalid_unsigned_int, const bool remove_duplicates=true, const double tol=1e-12, const bool ignore_unsupported_element_types=true)
Definition grid_in.cc:3402
void read_abaqus(std::istream &in, const bool apply_all_indicators_to_manifolds=false)
Definition grid_in.cc:1168
static std::string default_suffix(const Format format)
Definition grid_in.cc:4241
void read_xda(std::istream &in)
Definition grid_in.cc:1385
void read_comsol_mphtxt(std::istream &in)
Definition grid_in.cc:1449
static void skip_comment_lines(std::istream &in, const char comment_start)
Definition grid_in.cc:3967
void attach_triangulation(Triangulation< dim, spacedim > &tria)
Definition grid_in.cc:153
void read_msh(std::istream &in)
Definition grid_in.cc:2005
static Format parse_format(const std::string &format_name)
Definition grid_in.cc:4276
void read_vtu(std::istream &in)
Definition grid_in.cc:587
void read_tecplot(std::istream &in)
Definition grid_in.cc:3393
Format default_format
Definition grid_in.h:950
ExodusIIData read_exodusii(const std::string &filename, const bool apply_all_indicators_to_manifolds=false)
Definition grid_in.cc:3772
void read_dbmesh(std::istream &in)
Definition grid_in.cc:1224
void read_ucd(std::istream &in, const bool apply_all_indicators_to_manifolds=false)
Definition grid_in.cc:914
void read(std::istream &in, Format format=Default)
Definition grid_in.cc:4174
static void debug_output_grid(const std::vector< CellData< dim > > &cells, const std::vector< Point< spacedim > > &vertices, std::ostream &out)
Definition grid_in.cc:3993
static std::string get_format_names()
Definition grid_in.cc:4333
void read_unv(std::istream &in)
Definition grid_in.cc:614
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:2928
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
Point< 3 > vertices[4]
#define DeclException0(Exception0)
Definition exceptions.h:471
static ::ExceptionBase & ExcInvalidDBMeshFormat()
static ::ExceptionBase & ExcUnknownIdentifier(std::string arg1)
static ::ExceptionBase & ExcInvalidGMSHInput(std::string arg1)
static ::ExceptionBase & ExcInvalidDBMESHInput(std::string arg1)
static ::ExceptionBase & ExcDBMESHWrongDimension(int arg1)
static ::ExceptionBase & ExcInvalidVertexIndex(int arg1, int arg2)
static ::ExceptionBase & ExcNoTriangulationSelected()
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
static ::ExceptionBase & ExcUnknownElementType(int arg1)
static ::ExceptionBase & ExcGmshNoCellInformation(unsigned int arg1, unsigned int arg2)
static ::ExceptionBase & ExcInvalidVertexIndexGmsh(int arg1, int arg2, int arg3)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
static ::ExceptionBase & ExcUnknownSectionType(int arg1)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition exceptions.h:562
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516
static ::ExceptionBase & ExcGmshUnsupportedGeometry(int arg1)
static const unsigned int invalid_unsigned_int
Definition types.h:220
std::vector< std::vector< int > > id_to_sideset_ids
Definition grid_in.h:687