Reference documentation for deal.II version GIT relicensing-216-gcda843ca70 2024-03-28 12:20:02+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
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Private Member Functions | Private Attributes | List of all members
GridIn< dim, spacedim > Class Template Reference

#include <deal.II/grid/grid_in.h>

Classes

struct  ExodusIIData
 

Public Types

enum  Format {
  Default , unv , ucd , abaqus ,
  dbmesh , xda , msh , tecplot ,
  vtk , vtu , assimp , exodusii
}
 

Public Member Functions

 GridIn ()
 
 GridIn (Triangulation< dim, spacedim > &tria)
 
void attach_triangulation (Triangulation< dim, spacedim > &tria)
 
void read (std::istream &in, Format format=Default)
 
void read (const std::string &in, Format format=Default)
 
void read_vtk (std::istream &in)
 
void read_vtu (std::istream &in)
 
void read_unv (std::istream &in)
 
void read_ucd (std::istream &in, const bool apply_all_indicators_to_manifolds=false)
 
void read_abaqus (std::istream &in, const bool apply_all_indicators_to_manifolds=false)
 
void read_dbmesh (std::istream &in)
 
void read_xda (std::istream &in)
 
void read_msh (std::istream &in)
 
void read_msh (const std::string &filename)
 
void read_comsol_mphtxt (std::istream &in)
 
void read_tecplot (std::istream &in)
 
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)
 
ExodusIIData read_exodusii (const std::string &filename, const bool apply_all_indicators_to_manifolds=false)
 
void read_tecplot (std::istream &in)
 

Static Public Member Functions

static std::string default_suffix (const Format format)
 
static Format parse_format (const std::string &format_name)
 
static std::string get_format_names ()
 
static ::ExceptionBaseExcUnknownSectionType (int arg1)
 
static ::ExceptionBaseExcUnknownElementType (int arg1)
 
static ::ExceptionBaseExcUnknownIdentifier (std::string arg1)
 
static ::ExceptionBaseExcNoTriangulationSelected ()
 
static ::ExceptionBaseExcInvalidVertexIndex (int arg1, int arg2)
 
static ::ExceptionBaseExcInvalidVertexIndexGmsh (int arg1, int arg2, int arg3)
 
static ::ExceptionBaseExcInvalidDBMeshFormat ()
 
static ::ExceptionBaseExcInvalidDBMESHInput (std::string arg1)
 
static ::ExceptionBaseExcDBMESHWrongDimension (int arg1)
 
static ::ExceptionBaseExcInvalidGMSHInput (std::string arg1)
 
static ::ExceptionBaseExcGmshUnsupportedGeometry (int arg1)
 
static ::ExceptionBaseExcGmshNoCellInformation (unsigned int arg1, unsigned int arg2)
 

Protected Member Functions

void debug_output_grid (const std::vector< CellData< 2 > > &cells, const std::vector< Point< 2 > > &vertices, std::ostream &out)
 
void debug_output_grid (const std::vector< CellData< 3 > > &cells, const std::vector< Point< 3 > > &vertices, std::ostream &out)
 

Static Protected Member Functions

static void debug_output_grid (const std::vector< CellData< dim > > &cells, const std::vector< Point< spacedim > > &vertices, std::ostream &out)
 

Protected Attributes

SmartPointer< Triangulation< dim, spacedim >, GridIn< dim, spacedim > > tria
 

Static Private Member Functions

static void skip_empty_lines (std::istream &in)
 
static void skip_comment_lines (std::istream &in, const char comment_start)
 
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)
 

Private Attributes

Format default_format
 

Detailed Description

template<int dim, int spacedim = dim>
class GridIn< dim, spacedim >

This class implements an input mechanism for grid data. It allows to read a grid structure into a triangulation object. At present, UCD (unstructured cell data), DB Mesh, XDA, Gmsh, Tecplot, UNV, VTK, ASSIMP, and Cubit are supported as input format for grid data. Any numerical data other than geometric (vertex locations) and topological (how vertices form cells, faces, and edges) information is ignored, but the readers for the various formats generally do read information that associates material ids or boundary ids to cells or faces (see this and this glossary entry for more information).

Note
Since deal.II only supports line, quadrilateral and hexahedral meshes, the functions in this class can only read meshes that consist exclusively of such cells. If you absolutely need to work with a mesh that uses triangles or tetrahedra, then your only option is to convert the mesh to quadrilaterals and hexahedra. A tool that can do this is tethex, available here.

The mesh you read will form the coarsest level of a Triangulation object. As such, it must not contain hanging nodes or other forms of adaptive refinement, or strange things will happen if the mesh represented by the input file does in fact have them. This is due to the fact that most mesh description formats do not store neighborship information between cells, so the grid reading functions have to regenerate it. They do so by checking whether two cells have a common face. If there are hanging nodes in a triangulation, adjacent cells have no common (complete) face, so the grid reader concludes that the adjacent cells have no neighbors along these faces and must therefore be at the boundary. In effect, an internal crack of the domain is introduced this way. Since such cases are very hard to detect (how is GridIn supposed to decide whether a place where the faces of two small cells coincide with the face or a larger cell is in fact a hanging node associated with local refinement, or is indeed meant to be a crack in the domain?), the library does not make any attempt to catch such situations, and you will get a triangulation that probably does not do what you want. If your goal is to save and later read again a triangulation that has been adaptively refined, then this class is not your solution; rather take a look at the PersistentTriangulation class.

To read grid data, the triangulation to be filled has to be empty. Upon calling the functions of this class, the input file may contain only lines in one dimension; lines and quads in two dimensions; and lines, quads, and hexes in three dimensions. All other cell types (e.g. triangles in two dimensions, triangles or tetrahedra in 3d) are rejected. (Here, the "dimension" refers to the dimensionality of the mesh; it may be embedded in a higher dimensional space, such as a mesh on the two-dimensional surface of the sphere embedded in 3d, or a 1d mesh that discretizes a line in 3d.) The result will be a triangulation that consists of the cells described in the input file, and to the degree possible with material indicators and boundary indicators correctly set as described in the input file.

Note
You can not expect vertex and cell numbers in the triangulation to match those in the input file. (This is already clear based on the fact that we number cells and vertices separately, whereas this is not the case for some input file formats; some formats also do not require consecutive numbering, or start numbering at indices other than zero.)

Supported input formats

At present, the following input formats are supported:

Structure of input grid data.

It is your duty to use a correct numbering of vertices in the cell list, i.e. for lines in 1d, you have to first give the vertex with the lower coordinate value, then that with the higher coordinate value. For quadrilaterals in two dimensions, the vertex indices in the quad list have to be such that the vertices are numbered in counter-clockwise sense.

In two dimensions, another difficulty occurs, which has to do with the sense of a quadrilateral. A quad consists of four lines which have a direction, which is by definition as follows:

*   3-->--2
*   |     |
*   ^     ^
*   |     |
*   0-->--1
* 

Now, two adjacent cells must have a vertex numbering such that the direction of the common side is the same. For example, the following two quads

*   3---4---5
*   |   |   |
*   0---1---2
* 

may be characterised by the vertex numbers (0 1 4 3) and (1 2 5 4), since the middle line would get the direction 1->4 when viewed from both cells. The numbering (0 1 4 3) and (5 4 1 2) would not be allowed, since the left quad would give the common line the direction 1->4, while the right one would want to use 4->1, leading to an ambiguity. The Triangulation object is capable of detecting this special case, which can be eliminated by rotating the indices of the right quad by two. However, it would not know what to do if you gave the vertex indices (4 1 2 5), since then it would have to rotate by one element or three, the decision which to take is not yet implemented.

There are more ambiguous cases, where the triangulation may not know what to do at all without the use of sophisticated algorithms. Furthermore, similar problems exist in three space dimensions, where faces and lines have orientations that need to be taken care of.

For this reason, the read_* functions of this class that read in grids in various input formats call the GridTools::consistently_order_cells() function to bring the order of vertices that define the cells into an ordering that satisfies the requirements of the Triangulation class. Be sure to read the documentation of that class if you experience unexpected problems when reading grids through this class.

Dealing with distorted mesh cells

For each of the mesh reading functions, the last call is always to Triangulation::create_triangulation(). That function checks whether all the cells it creates as part of the coarse mesh are distorted or not (where distortion here means that the Jacobian of the mapping from the reference cell to the real cell has a non-positive determinant, i.e. the cell is pinched or twisted; see the entry on distorted cells in the glossary). If it finds any such cells, it throws an exception. This exception is not caught in the grid reader functions of the current class, and so will propagate through to the function that called it. There, you can catch and ignore the exception if you are certain that there is no harm in dealing with such cells. If you were not aware that your mesh had such cells, your results will likely be of dubious quality at best if you ignore the exception.

Definition at line 310 of file grid_in.h.

Member Enumeration Documentation

◆ Format

template<int dim, int spacedim = dim>
enum GridIn::Format

List of possible mesh input formats. These values are used when calling the function read() in order to determine the actual reader to be called.

Enumerator
Default 

Use GridIn::default_format stored in this object.

unv 

Use read_unv()

ucd 

Use read_ucd()

abaqus 

Use read_abaqus()

dbmesh 

Use read_dbmesh()

xda 

Use read_xda()

msh 

Use read_msh()

tecplot 

Use read_tecplot()

vtk 

Use read_vtk()

vtu 

Use read_vtu()

assimp 

Use read_assimp()

exodusii 

Use read_exodusii()

Definition at line 317 of file grid_in.h.

Constructor & Destructor Documentation

◆ GridIn() [1/2]

template<int dim, int spacedim>
GridIn< dim, spacedim >::GridIn ( )

Constructor.

Definition at line 137 of file grid_in.cc.

◆ GridIn() [2/2]

template<int dim, int spacedim>
GridIn< dim, spacedim >::GridIn ( Triangulation< dim, spacedim > &  tria)

Constructor. Attach this triangulation to be fed with the grid data.

Definition at line 145 of file grid_in.cc.

Member Function Documentation

◆ attach_triangulation()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::attach_triangulation ( Triangulation< dim, spacedim > &  tria)

Attach this triangulation to be fed with the grid data.

Definition at line 154 of file grid_in.cc.

◆ read() [1/2]

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read ( std::istream &  in,
Format  format = Default 
)

Read from the given stream. If no format is given, GridIn::Format::Default is used.

Definition at line 4182 of file grid_in.cc.

◆ read() [2/2]

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read ( const std::string &  in,
Format  format = Default 
)

Open the file given by the string and call the previous function read(). This function uses the PathSearch mechanism to find files. The file class used is MESH.

Definition at line 4140 of file grid_in.cc.

◆ read_vtk()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_vtk ( std::istream &  in)

Read grid data from a unstructured vtk file. The vtk file may contain the following VTK cell types: VTK_HEXAHEDRON (12), VTK_TETRA (10), VTK_QUAD (9), VTK_TRIANGLE (5), and VTK_LINE (3).

Depending on the template dimension, only some of the above are accepted.

In particular, in three dimensions, this function expects the file to contain

  • VTK_HEXAHEDRON/VTK_TETRA cell types
  • VTK_QUAD/VTK_TRIANGLE cell types, to specify optional boundary or interior quad faces
  • VTK_LINE cell types, to specify optional boundary or interior edges

In two dimensions:

  • VTK_QUAD/VTK_TRIANGLE cell types
  • VTK_LINE cell types, to specify optional boundary or interior edges

In one dimension

  • VTK_LINE cell types

The input file may specify boundary ids, material ids, and manifold ids using the CELL_DATA section of the VTK file format.

This function interprets two types of CELL_DATA contained in the input file: SCALARS MaterialID, used to specify the material_id of the cells, or the boundary_id of the faces and edges, and SCALARS ManifoldID, that can be used to specify the manifold id of any Triangulation object (cell, face, or edge).

The companion GridOut::write_vtk function can be used to write VTK files compatible with this method.

Also see Simplex support.

Definition at line 163 of file grid_in.cc.

◆ read_vtu()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_vtu ( std::istream &  in)

Read grid data from a unstructured vtu file, saved by deal.II using GridOut::write_vtu(), with the flag GridOutFlags::Vtu::serialize_triangulation set to true.

Notice that this function does not support reading in arbitrary vtu files, but only files that were written by deal.II itself, using the function GridOut::write_vtu and setting GridOutFlags::Vtu::serialize_triangulation to true.

When this flag is set to true, the generated vtu file contains the triangulation in a xml section which is ignored by general vtu readers. If this section is absent, an exception is thrown.

Definition at line 588 of file grid_in.cc.

◆ read_unv()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_unv ( std::istream &  in)

Read grid data from an unv file as generated by the Salome mesh generator. Numerical data is ignored.

Note the comments on generating this file format in the general documentation of this class.

Definition at line 615 of file grid_in.cc.

◆ read_ucd()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_ucd ( std::istream &  in,
const bool  apply_all_indicators_to_manifolds = false 
)

Read grid data from an ucd file. Numerical data is ignored. It is not possible to use a ucd file to set both boundary_id and manifold_id for the same cell. Yet it is possible to use the flag apply_all_indicators_to_manifolds to decide if the indicators in the file refer to manifolds (flag set to true) or boundaries (flag set to false). If the flag is set, the indicators are used for cells as manifold id, too.

Definition at line 915 of file grid_in.cc.

◆ read_abaqus()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_abaqus ( std::istream &  in,
const bool  apply_all_indicators_to_manifolds = false 
)

Read grid data from an Abaqus file. Numerical and constitutive data is ignored. As in the case of the ucd file format, it is possible to use the flag apply_all_indicators_to_manifolds to decide if the indicators in the file refer to manifolds (flag set to true) or boundaries (flag set to false).

Note
The current implementation of this mesh reader is suboptimal, and may therefore be slow for large meshes.
Usage tips for Cubit:
  • Multiple material-id's can be defined in the mesh. This is done by specifying blocksets in the pre-processor.
  • Arbitrary surface boundaries can be defined in the mesh. This is done by specifying sidesets in the pre-processor. In particular, boundaries are not confined to just surfaces (in 3d) individual element faces can be added to the sideset as well. This is useful when a boundary condition is to be applied on a complex shape boundary that is difficult to define using "surfaces" alone. Similar can be done in 2d.
Compatibility information for this file format is listed below.
  • Files generated in Abaqus CAE 6.12 have been verified to be correctly imported, but older (or newer) versions of Abaqus may also generate valid input decks.
  • Files generated using Cubit 11.x, 12.x, 13.x, 14.x and 15.x are valid, but only when using a specific set of export steps. These are as follows:
    • Go to "Analysis setup mode" by clicking on the disc icon in the toolbar on the right.
    • Select "Export Mesh" under "Operation" by clicking on the necessary icon in the toolbar on the right.
    • Select an output file. In Cubit version 11.0 and 12.0 it might be necessary to click on the browse button and type it in the dialogue that pops up.
    • Select the dimension to output in.
    • Tick the overwrite box.
    • If using Cubit v12.0 onwards, uncheck the box "Export using Cubit ID's". An invalid file will encounter errors if this box is left checked.
    • Click apply.

Definition at line 1169 of file grid_in.cc.

◆ read_dbmesh()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_dbmesh ( std::istream &  in)

Read grid data from a file containing data in the DB mesh format.

Definition at line 1225 of file grid_in.cc.

◆ read_xda()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_xda ( std::istream &  in)

Read grid data from a file containing data in the XDA format.

Definition at line 1386 of file grid_in.cc.

◆ read_msh() [1/2]

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_msh ( std::istream &  in)

Read grid data from an msh file. The Gmsh formats are documented at http://www.gmsh.info/.

Also see Simplex support.

Definition at line 2006 of file grid_in.cc.

◆ read_msh() [2/2]

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_msh ( const std::string &  filename)

Read grid data using Gmsh API. Any file supported by Gmsh can be passed as argument. The format is deduced from the filename extension.

This function interprets non-named physical ids (gmsh format < 4.0) as material or boundary ids (similarly to what happens with the other read_msh() function). If you want to specify non default manifold or boundary ids, you must group all entities that require a non default boundary or manifold id into named physical groups, where the name is interpreted using the function Patterns::Tools::to_value() applied to a std::map<std::string, int>. The keys can be either MaterialID (if the physical group refers to object of dimension dim), BoundaryID (if the group refers to objects of dimension < dim), or ManifoldID.

From the Gmsh documentation, the formats of the physical tags follows the following conventions:

\$PhysicalNames // same as MSH version 2
numPhysicalNames(ASCII int)
dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max)
...
\$EndPhysicalNames

For example, the following snippet of mesh file

MeshFormat
4.1 0 8
\$EndMeshFormat
\$PhysicalNames
4
1 1 "ManifoldID:0"
1 2 "BoundaryID: -1, ManifoldID: 1"
2 3 "ManifoldID: 1"
2 4 "MaterialID: 2, ManifoldID: 1"
\$EndPhysicalNames
\$Entities
...

refers to a two dimensional grid where:

  • a portion of the boundary of dimension 1 has physical tag 1, and manifold id 0
  • some internal faces (lines of dimension 1) have manifold id 1
  • some elements have manifold id 1 (and material id equal to the default value, i.e., zero)
  • some elements have manifold id 1 and material id equal to 2

If the physical groups are not named, then the behavior is the same as the other read_msh() function, i.e., the physical tag itself is interpreted as a boundary or material id. Physical surface numbers created in Gmsh, which can be seen in the .geo file, become material IDs.

Also see Simplex support.

Definition at line 2695 of file grid_in.cc.

◆ read_comsol_mphtxt()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_comsol_mphtxt ( std::istream &  in)

Read grid data from a .mphtxt file. .mphtxt is one of the file formats typically generated by COMSOL. The file format is described at http://victorsndvg.github.io/FEconv/formats/mphtxt.xhtml .

The reader interprets the "geometric entity indicators" that COMSOL writes into these files as either boundary indicators (for edges and faces of cells) or as material ids (for cells). See the glossary for a description of boundary indicators. and material indicators

COMSOL has a habit of assigning "geometric entity indicators" not only to edges and faces on the actual boundary, but also to interior faces and edges. For example, for the following volume mesh generated by COMSOL,

the marked edges and faces are as follows:

Here, some of the marked lines and faces with explicitly given geometric entity indicators are in the interior of the domain – an artifact of the geometry description that was used to describe the mesh. However, we can of course not assign boundary indicators to interior edges and faces. As a consequence, this reader function simply ignores the geometric entity indicator for edges and faces that are not in fact on the boundary of the domain. The result is then a mesh in which only the following edges and faces are explicitly assigned boundary indicators:

Also see Simplex support.

Definition at line 1450 of file grid_in.cc.

◆ read_tecplot() [1/2]

template<int dim, int spacedim>
void GridIn< dim, spacedim >::read_tecplot ( std::istream &  in)

Read grid data from a file containing tecplot ASCII data. This also works in the absence of any tecplot installation.

Definition at line 3394 of file grid_in.cc.

◆ read_assimp()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::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 
)

Read in a file supported by Assimp, and generate a Triangulation out of it. If you specify a mesh_index, only the mesh with the given index will be extracted, otherwise all meshes which are present in the file will be used to generate the Triangulation.

This function can only be used to read two-dimensional meshes (possibly embedded in three dimensions). This is the standard for graphical software such as blender, or 3d studio max, and that is what the original Assimp library was built for. We "bend" it to deal.II to support complex co-dimension one meshes and complex two-dimensional meshes.

If remove_duplicates is set to true (the default), then duplicated vertices will be removed if their distance is lower than tol.

Only the elements compatible with the given dimension and space dimension will be extracted from the mesh, and only those elements that are compatible with deal.II are supported. If you set ignore_unsupported_element_types, all the other element types are simply ignored by this algorithm. If your mesh contains a mixture of triangles and quadrilaterals, for example, only the quadrilaterals will be extracted. The resulting mesh (as represented in the Triangulation object) may not make any sense if you are mixing compatible and incompatible element types. If ignore_unsupported_element_types is set to false, then an exception is thrown when an unsupported type is encountered.

Parameters
filenameThe file to read from
mesh_indexIndex of the mesh within the file
remove_duplicatesRemove duplicated vertices
tolTolerance to use when removing vertices
ignore_unsupported_element_typesDon't throw exceptions if we encounter unsupported types during parsing

Definition at line 3403 of file grid_in.cc.

◆ read_exodusii()

template<int dim, int spacedim>
GridIn< dim, spacedim >::ExodusIIData GridIn< dim, spacedim >::read_exodusii ( const std::string &  filename,
const bool  apply_all_indicators_to_manifolds = false 
)

Read in a mesh stored in the ExodusII file format.

ExodusII is a feature-rich file format that supports many more features (like node sets, finite element fields, quality assurance data, and more) than most other grid formats supported by this class. Many of these features do not have equivalent representations in deal.II and are therefore not supported (for example, deal.II does not assign degrees of freedom directly to nodes, so data stored in a nodal format is not loaded by this function). At the current time only the following information is extracted from the input file :

  1. Block ids: the block id of an element is loaded as its material id.
  2. Elements and vertices: the core geometric information stored in the ExodusII file populates the attached Triangulation object. Higher-order elements are automatically truncated to lower-order elements since deal.II does not support this feature (e.g., there is no equivalent to the QUAD9 element in deal.II since all quadrilaterals have four vertices and additional geometric information is either stored in a Manifold or something like MappingQEulerian).
  3. Sideset ids: these are interpreted as boundary ids or manifold ids (see the note on the output value below). An error will occur if you attempt to read an ExodusII file that assigns a sideset id to an internal face boundary id.

Sideset ids are not translated for Triangulations with nonzero codimension since those Triangulations do not support the setting of boundary ids.

Parameters
filenameThe name of the file to read from.
apply_all_indicators_to_manifoldsBoolean determining if the sideset ids should be interpreted as manifold ids or boundary ids. The default value is false, i.e., treat all sideset ids as boundary ids. If your mesh sets sideset ids on internal faces then it will be necessary to set this argument to true and then do some postprocessing to set the boundary ids correctly.
Returns
This function returns a struct containing some extra data stored by the ExodusII file that cannot be loaded into a Triangulation - see ExodusIIData for more information.

A cell face in ExodusII can be in an arbitrary number of sidesets (i.e., it can have an arbitrary number of sideset ids) - however, a boundary cell face in deal.II has exactly one boundary id. All boundary faces that are not in a sideset are given the (default) boundary id of \(0\). This function then groups sidesets together into unique sets and gives each one a boundary id. For example: Consider a single-quadrilateral mesh whose left side has no sideset id, right side has sideset ids \(0\) and \(1\), and whose bottom and top sides have sideset ids of \(0\). The left face will have a boundary id of \(0\), the top and bottom faces boundary ids of \(1\), and the right face a boundary id of \(2\). Hence the vector returned by this function in that case will be \(\{\{\}, \{0\}, \{0, 1\}\}\).

Definition at line 3773 of file grid_in.cc.

◆ default_suffix()

template<int dim, int spacedim>
std::string GridIn< dim, spacedim >::default_suffix ( const Format  format)
static

Return the standard suffix for a file in this format.

Definition at line 4249 of file grid_in.cc.

◆ parse_format()

template<int dim, int spacedim>
GridIn< dim, spacedim >::Format GridIn< dim, spacedim >::parse_format ( const std::string &  format_name)
static

Return the enum Format for the format name.

Definition at line 4284 of file grid_in.cc.

◆ get_format_names()

template<int dim, int spacedim>
std::string GridIn< dim, spacedim >::get_format_names ( )
static

Return a list of implemented input formats. The different names are separated by vertical bar signs (‘|’) as used by the ParameterHandler classes.

Definition at line 4341 of file grid_in.cc.

◆ debug_output_grid() [1/3]

template<int dim, int spacedim>
void GridIn< dim, spacedim >::debug_output_grid ( const std::vector< CellData< dim > > &  cells,
const std::vector< Point< spacedim > > &  vertices,
std::ostream &  out 
)
staticprotected

This function can write the raw cell data objects created by the read_* functions in Gnuplot format to a stream. This is sometimes handy if one would like to see what actually was created, if it is known that the data is not correct in some way, but the Triangulation class refuses to generate a triangulation because of these errors. In particular, the output of this class writes out the cell numbers along with the direction of the faces of each cell. In particular the latter information is needed to verify whether the cell data objects follow the requirements of the ordering of cells and their faces, i.e. that all faces need to have unique directions and specified orientations with respect to neighboring cells (see the documentations to this class and the GridTools::consistently_order_cells() function).

The output of this function consists of vectors for each line bounding the cells indicating the direction it has with respect to the orientation of this cell, and the cell number. The whole output is in a form such that it can be read in by Gnuplot and generate the full plot without further ado by the user.

Definition at line 3994 of file grid_in.cc.

◆ skip_empty_lines()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::skip_empty_lines ( std::istream &  in)
staticprivate

Skip empty lines in the input stream, i.e. lines that contain either nothing or only whitespace.

Definition at line 3938 of file grid_in.cc.

◆ skip_comment_lines()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::skip_comment_lines ( std::istream &  in,
const char  comment_start 
)
staticprivate

Skip lines of comment that start with the indicated character (e.g. #) following the point where the given input stream presently is. After the call to this function, the stream is at the start of the first line after the comment lines, or at the same position as before if there were no lines of comments.

Definition at line 3968 of file grid_in.cc.

◆ parse_tecplot_header()

template<int dim, int spacedim>
void GridIn< dim, spacedim >::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 
)
staticprivate

This function does the nasty work (due to very lax conventions and different versions of the tecplot format) of extracting the important parameters from a tecplot header, contained in the string header. The other variables are output variables, their value has no influence on the function execution..

Definition at line 2929 of file grid_in.cc.

◆ read_tecplot() [2/2]

void GridIn< 2 >::read_tecplot ( std::istream &  in)

Definition at line 3163 of file grid_in.cc.

◆ debug_output_grid() [2/3]

void GridIn< 2 >::debug_output_grid ( const std::vector< CellData< 2 > > &  cells,
const std::vector< Point< 2 > > &  vertices,
std::ostream &  out 
)
protected

Definition at line 4006 of file grid_in.cc.

◆ debug_output_grid() [3/3]

void GridIn< 3 >::debug_output_grid ( const std::vector< CellData< 3 > > &  cells,
const std::vector< Point< 3 > > &  vertices,
std::ostream &  out 
)
protected

Definition at line 4067 of file grid_in.cc.

Member Data Documentation

◆ tria

template<int dim, int spacedim = dim>
SmartPointer<Triangulation<dim, spacedim>, GridIn<dim, spacedim> > GridIn< dim, spacedim >::tria
protected

Store address of the triangulation to be fed with the data read in.

Definition at line 884 of file grid_in.h.

◆ default_format

template<int dim, int spacedim = dim>
Format GridIn< dim, spacedim >::default_format
private

Input format used by read() if no format is given.

Definition at line 949 of file grid_in.h.


The documentation for this class was generated from the following files: