Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08: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
Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members

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

Public Types

enum  OutputFormat {
  none , dx , gnuplot , eps ,
  ucd , xfig , msh , svg ,
  mathgl , vtk , vtu
}
 

Public Member Functions

 GridOut ()
 
template<int dim, int spacedim>
void write_dx (const Triangulation< dim, spacedim > &tria, std::ostream &out) const
 
template<int dim, int spacedim>
void write_gnuplot (const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
 
template<int dim, int spacedim>
void write_msh (const Triangulation< dim, spacedim > &tria, std::ostream &out) const
 
template<int dim, int spacedim>
void write_msh (const Triangulation< dim, spacedim > &tria, const std::string &filename) const
 
template<int dim, int spacedim>
void write_ucd (const Triangulation< dim, spacedim > &tria, std::ostream &out) const
 
template<int dim, int spacedim>
void write_eps (const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
 
template<int dim, int spacedim>
void write_xfig (const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
 
void write_svg (const Triangulation< 2, 2 > &tria, std::ostream &out) const
 
template<int dim, int spacedim>
void write_svg (const Triangulation< dim, spacedim > &tria, std::ostream &out) const
 
template<int dim, int spacedim>
void write_mathgl (const Triangulation< dim, spacedim > &tria, std::ostream &out) const
 
template<int dim, int spacedim>
void write_vtk (const Triangulation< dim, spacedim > &tria, std::ostream &out) const
 
template<int dim, int spacedim>
void write_vtu (const Triangulation< dim, spacedim > &tria, std::ostream &out) const
 
template<int dim, int spacedim>
void write_mesh_per_processor_as_vtu (const Triangulation< dim, spacedim > &tria, const std::string &filename_without_extension, const bool view_levels=false, const bool include_artificial=false) const
 
template<int dim, int spacedim>
void write (const Triangulation< dim, spacedim > &tria, std::ostream &out, const OutputFormat output_format, const Mapping< dim, spacedim > *mapping=nullptr) const
 
template<int dim, int spacedim>
void write (const Triangulation< dim, spacedim > &tria, std::ostream &out, const Mapping< dim, spacedim > *mapping=nullptr) const
 
void set_flags (const GridOutFlags::DX &flags)
 
void set_flags (const GridOutFlags::Msh &flags)
 
void set_flags (const GridOutFlags::Ucd &flags)
 
void set_flags (const GridOutFlags::Gnuplot &flags)
 
void set_flags (const GridOutFlags::Eps< 1 > &flags)
 
void set_flags (const GridOutFlags::Eps< 2 > &flags)
 
void set_flags (const GridOutFlags::Eps< 3 > &flags)
 
void set_flags (const GridOutFlags::XFig &flags)
 
void set_flags (const GridOutFlags::Svg &flags)
 
void set_flags (const GridOutFlags::MathGL &flags)
 
void set_flags (const GridOutFlags::Vtk &flags)
 
void set_flags (const GridOutFlags::Vtu &flags)
 
std::string default_suffix () const
 
void parse_parameters (ParameterHandler &param)
 
std::size_t memory_consumption () const
 
template<>
void write_dx (const Triangulation< 1 > &, std::ostream &) const
 
template<>
void write_dx (const Triangulation< 1, 2 > &, std::ostream &) const
 
template<>
void write_dx (const Triangulation< 1, 3 > &, std::ostream &) const
 
template<>
void write_xfig (const Triangulation< 2 > &tria, std::ostream &out, const Mapping< 2 > *) const
 
template<>
void write_mathgl (const Triangulation< 1 > &, std::ostream &) const
 

Static Public Member Functions

static std::string default_suffix (const OutputFormat output_format)
 
static OutputFormat parse_output_format (const std::string &format_name)
 
static std::string get_output_format_names ()
 
static void declare_parameters (ParameterHandler &param)
 
static ::ExceptionBaseExcInvalidState ()
 

Private Member Functions

template<int dim, int spacedim>
unsigned int write_msh_faces (const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_msh_faces (const Triangulation< 1, 1 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_msh_faces (const Triangulation< 1, 2 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_msh_faces (const Triangulation< 1, 3 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
template<int dim, int spacedim>
unsigned int write_msh_lines (const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_msh_lines (const Triangulation< 1, 1 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_msh_lines (const Triangulation< 1, 2 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_msh_lines (const Triangulation< 1, 3 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_msh_lines (const Triangulation< 2, 2 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_msh_lines (const Triangulation< 2, 3 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
template<int dim, int spacedim>
unsigned int write_ucd_faces (const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_ucd_faces (const Triangulation< 1, 1 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_ucd_faces (const Triangulation< 1, 2 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_ucd_faces (const Triangulation< 1, 3 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
template<int dim, int spacedim>
unsigned int write_ucd_lines (const Triangulation< dim, spacedim > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_ucd_lines (const Triangulation< 1, 1 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_ucd_lines (const Triangulation< 1, 2 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_ucd_lines (const Triangulation< 1, 3 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_ucd_lines (const Triangulation< 2, 2 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
unsigned int write_ucd_lines (const Triangulation< 2, 3 > &tria, const unsigned int next_element_index, std::ostream &out) const
 
template<int dim, int spacedim>
unsigned int n_boundary_faces (const Triangulation< dim, spacedim > &tria) const
 
unsigned int n_boundary_faces (const Triangulation< 1, 1 > &tria) const
 
unsigned int n_boundary_faces (const Triangulation< 1, 2 > &tria) const
 
unsigned int n_boundary_faces (const Triangulation< 1, 3 > &tria) const
 
template<int dim, int spacedim>
unsigned int n_boundary_lines (const Triangulation< dim, spacedim > &tria) const
 
unsigned int n_boundary_lines (const Triangulation< 1, 1 > &tria) const
 
unsigned int n_boundary_lines (const Triangulation< 1, 2 > &tria) const
 
unsigned int n_boundary_lines (const Triangulation< 1, 3 > &tria) const
 
unsigned int n_boundary_lines (const Triangulation< 2, 2 > &tria) const
 
unsigned int n_boundary_lines (const Triangulation< 2, 3 > &tria) const
 

Private Attributes

OutputFormat default_format
 
GridOutFlags::DX dx_flags
 
GridOutFlags::Msh msh_flags
 
GridOutFlags::Ucd ucd_flags
 
GridOutFlags::Gnuplot gnuplot_flags
 
GridOutFlags::Eps< 1 > eps_flags_1
 
GridOutFlags::Eps< 2 > eps_flags_2
 
GridOutFlags::Eps< 3 > eps_flags_3
 
GridOutFlags::XFig xfig_flags
 
GridOutFlags::Svg svg_flags
 
GridOutFlags::MathGL mathgl_flags
 
GridOutFlags::Vtk vtk_flags
 
GridOutFlags::Vtu vtu_flags
 

Detailed Description

This class provides a means to output a triangulation to a file in different formats. See the enum GridOut::OutputFormat for a list of formats and the corresponding output function names.

Usage is simple: either you use the direct form

std::ofstream output_file("some_filename");
GridOut().write_gnuplot (tria, output_file);

if you know which format you want to have, or if you want the format to be a runtime parameter, you can write

GridOut::OutputFormat output_format =
GridOut::parse_output_format(get_format_name_from_somewhere());
std::ofstream output_file("some_filename"
+ GridOut::default_suffix(output_format));
GridOut().write (tria, output_file, output_format);
std::string default_suffix() const
Definition grid_out.cc:593
static OutputFormat parse_output_format(const std::string &format_name)
Definition grid_out.cc:601

The function get_output_format_names() provides a list of possible names of output formats in a string that is understandable by the ParameterHandler class.

Note that here, we have created an unnamed object of type GridOut and called one of its write_* functions. This looks like as if the respective function could really be made static. This was not done in order to allow for parameters to be passed to the different output functions in a way compatible with the scheme of allowing the right output format to be selected at run-time through the generic write function.

In order to explain this, consider each function had one or more additional parameters giving the details of output, for example position of the spectator for 3d meshed, line thicknesses, etc. While this would allow each output function any flexibility it needs, it would not allow us to use the generic function write which is given a parameter determining the output format, since it is impractical to give it a list of parameters for each and every output format supported which it may then pass on to the respective output function.

Rather, we have chosen to let each object of this class GridOut have a set of parameters for each supported output format. These are collected in structures GridOutFlags::Eps(), GridOutFlags::Gnuplot(), etc declared in the GridOutFlags namespace, and you can set your preferred flags like this:

GridOut grid_out;
... // set some fields in ucd_flags
grid_out.set_flags (ucd_flags);
...
... // write some file with data_out
void set_flags(const GridOutFlags::DX &flags)
Definition grid_out.cc:471
GridOutFlags::Ucd ucd_flags
Definition grid_out.h:1566

The respective output function then use the so-set flags. By default, they are set to reasonable values as described above and in the documentation of the different flags structures. Resetting the flags can be done by calling grid_out.set_flags (GridOutFlags::Ucd());, since the default constructor of each of the flags structures sets the parameters to their initial values.

The advantage of this approach is that it is possible to change the flags of one or more output formats according to your needs and later use the generic write function; the actual output function then called will use the flags as set before.

Note that some of the structures describing the flags of the different output formats are empty since the respective format does not support any flags. The structure and the set_flags function are provided anyway. Note also that some of the structures may differ between the dimensions supported by this class; they then have a template parameter, as usual.

Definition at line 992 of file grid_out.h.

Member Enumeration Documentation

◆ OutputFormat

Declaration of a name for each of the different output formats. These are used by the generic output function write() to determine the actual output format.

Enumerator
none 

Do nothing in write()

dx 

write() calls write_dx()

gnuplot 

write() calls write_gnuplot()

eps 

write() calls write_eps()

ucd 

write() calls write_ucd()

xfig 

write() calls write_xfig()

msh 

write() calls write_msh()

svg 

write() calls write_svg()

mathgl 

write() calls write_mathgl()

vtk 

write() calls write_vtk()

vtu 

write() calls write_vtu()

Definition at line 1000 of file grid_out.h.

Constructor & Destructor Documentation

◆ GridOut()

GridOut::GridOut ( )

Constructor.

Definition at line 465 of file grid_out.cc.

Member Function Documentation

◆ write_dx() [1/4]

template<int dim, int spacedim>
void GridOut::write_dx ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out 
) const

Write triangulation in OpenDX format.

Cells or faces are written together with their level and their material id or boundary indicator, resp.

Not implemented for the codimension one case.

Definition at line 781 of file grid_out.cc.

◆ write_gnuplot()

template<int dim, int spacedim>
void GridOut::write_gnuplot ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out,
const Mapping< dim, spacedim > *  mapping = nullptr 
) const

Write the triangulation in the gnuplot format.

In GNUPLOT format, each cell is written as a sequence of its confining lines. Apart from the coordinates of the lines' end points, the level and the material of the cell are appended to each line of output. Therefore, if you let GNUPLOT draw a 2d grid as a 3d plot, you will see more refined cells being raised against cells with less refinement. Also, if you draw a cut through a 3d grid, you can extrude the refinement level in the direction orthogonal to the cut plane. The same can be done with the material id, which is plotted after the level.

A more useful application of this feature is the following: if you use the GNUPLOT command (for a 2d grid here)

* splot [:][:][2.5:3.5] "grid_file.gnuplot" *
* 

then the whole x- and y-range will be plotted, i.e. the whole grid, but only those lines with a z-value between 2.5 and 3.5. Since the z-values were chosen to be the level to which a cell belongs, this results in a plot of those cells only that belong to level 3 in this example. This way, it is easy to produce plots of the different levels of grid.

mapping is a pointer to a mapping used for the transformation of cells at the boundary. If zero, then use standard Q1 mapping.

Names and values of additional flags controlling the output can be found in the documentation of the GridOutFlags::Gnuplot class, which also describes some caveats for the codimension one case.

Definition at line 4598 of file grid_out.cc.

◆ write_msh() [1/2]

template<int dim, int spacedim>
void GridOut::write_msh ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out 
) const

Write the triangulation in the msh format.

Msh is the format used by Gmsh and it is described in the Gmsh user's guide. Besides the usual output of the grid only, you can decide through additional flags (see below, and the documentation of the GridOutFlags::Msh() class) whether boundary faces with non-zero boundary indicator shall be written to the file explicitly. This is useful, if you want to re-read the grid later on, since deal.II sets the boundary indicator to zero by default; therefore, to obtain the same triangulation as before, you have to specify faces with differing boundary indicators explicitly, which is done by this flag.

Names and values of further flags controlling the output can be found in the documentation of the GridOutFlags::Msh() class.

Works also in the codimension one case.

Definition at line 1021 of file grid_out.cc.

◆ write_msh() [2/2]

template<int dim, int spacedim>
void GridOut::write_msh ( const Triangulation< dim, spacedim > &  tria,
const std::string &  filename 
) const

Write the triangulation in any format supported by gmsh API.

Gmsh API allows writing its output in several formats through their C++ API. This function translates a Triangulation object into a gmsh collection of entities and calls the gmsh::write() method with the filename passed as argument. This method generates a different entity for each unique pair of non default manifold id and boundary id, and writes a gmsh physical group for each unique combination, allowing you to read back the triangulation using the GridIn::read_msh() method that takes a string as argument.

In particular, all cell objects with non default boundary id or non default manifold id are grouped in a unique physical tag, whose name contains the boundary and manifold indicators. The names are constructed using Patterns::Tools::to_value() applied to a std::map<std::string, int> where the keys are either MaterialID, BoundaryID, or ManifoldID, i.e., a cell with material id 1, and manifold id 3 would be grouped in a physical tag (whose numbering is not specified), named MaterialID:1, ManifoldID:3.

For example, calling the method with a hyper ball grid refined once, would results in the following physical tags defined in the output file :

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

The special boundary id -1 is used to indicate internal boundaries. The internal boundaries must be specified whenever it is necessary to specify a non-flat manifold id.

Definition at line 1447 of file grid_out.cc.

◆ write_ucd()

template<int dim, int spacedim>
void GridOut::write_ucd ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out 
) const

Write the triangulation in the ucd format.

UCD (unstructured cell data) is the format used by AVS and some other programs. It is described in the AVS developer's guide. Besides the usual output of the grid only, you can decide through additional flags (see below, and the documentation of the GridOutFlags::Ucd() class) whether boundary faces with non-zero boundary indicator shall be written to the file explicitly. This is useful, if you want to re-read the grid later on, since deal.II sets the boundary indicator to zero by default; therefore, to obtain the same triangulation as before, you have to specify faces with differing boundary indicators explicitly, which is done by this flag.

Names and values of further flags controlling the output can be found in the documentation of the GridOutFlags::Ucd() class.

Works also for the codimension one case.

Definition at line 1126 of file grid_out.cc.

◆ write_eps()

template<int dim, int spacedim>
void GridOut::write_eps ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out,
const Mapping< dim, spacedim > *  mapping = nullptr 
) const

Write the triangulation in the encapsulated postscript format.

In this format, each line of the triangulation is written separately. We scale the picture such that either x-values or y-values range between zero and a fixed size. The other axis is scaled by the same factor. Which axis is taken to compute the scale and the size of the box it shall fit into is determined by the output flags (see below, and the documentation of the GridOutFlags::Eps() class).

The bounding box is close to the triangulation on all four sides, without an extra frame. The line width is chosen to be 0.5 by default, but can be changed. The line width is to be compared with the extension of the picture, of which the default is 300.

The flag color_lines_on_user_flag allows to draw lines with the user_flag set to be drawn in red. The colors black and red are defined as b and r in the preamble of the output file and can be changed there according to need.

mapping is a pointer to a mapping used for the transformation of cells at the boundary. If zero, then use standard Q1 mapping.

Names and values of additional flags controlling the output can be found in the documentation of the GridOutFlags::Eps() class. Especially the viewpoint for three dimensional grids is of importance here.

Not implemented for the codimension one case.

Definition at line 5109 of file grid_out.cc.

◆ write_xfig() [1/2]

template<int dim, int spacedim>
void GridOut::write_xfig ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out,
const Mapping< dim, spacedim > *  mapping = nullptr 
) const

Write two-dimensional XFig-file.

This function writes all grid cells as polygons and optionally boundary lines. Several parameters can be adjusted by the XFigFlags control object.

If levels are coded to depth, the complete grid hierarchy is plotted with fine cells before their parents. This way, levels can be switched on and off in xfig by selecting levels.

Polygons are either at depth 900-level or at 900+material_id, depending on the flag level_depth. Accordingly, boundary edges are at depth 800-level or at 800+boundary_id. Therefore, boundary edges are always in front of cells.

Not implemented for the codimension one case.

Definition at line 1243 of file grid_out.cc.

◆ write_svg() [1/2]

void GridOut::write_svg ( const Triangulation< 2, 2 > &  tria,
std::ostream &  out 
) const

Write the triangulation in the SVG format.

SVG (Scalable Vector Graphics) is an XML-based vector image format developed and maintained by the World Wide Web Consortium (W3C). This function conforms to the latest specification SVG 1.1, released on August 16, 2011.

The cells of the triangulation are written as polygons with additional lines at the boundary of the triangulation. A coloring of the cells is further possible in order to visualize a certain property of the cells such as their level or material id. A colorbar can be drawn to encode the chosen coloring. Moreover, a cell label can be added, showing level index, etc. Indeed, by using the set_flags() with an appropriately generated object of type GridOutFlags::Svg, many aspects of how and what is being visualized by this function can be customized.

Note
This function is currently only implemented for two-dimensional grids in two space dimensions.

Definition at line 1702 of file grid_out.cc.

◆ write_svg() [2/2]

template<int dim, int spacedim>
void GridOut::write_svg ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out 
) const

Declaration of the same function as above for all other dimensions and space dimensions. This function is not currently implemented and is only declared to exist to support dimension independent programming.

Definition at line 1684 of file grid_out.cc.

◆ write_mathgl() [1/2]

template<int dim, int spacedim>
void GridOut::write_mathgl ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out 
) const

Write triangulation in MathGL script format. To interpret this file a version of MathGL>=2.0.0 is required.

To get a handle on the resultant MathGL script within a graphical environment an interpreter is needed. A suggestion to start with is mglview, which is bundled with MathGL. mglview can interpret and display small-to-medium MathGL scripts in a graphical window and enables conversion to other formats such as EPS, PNG, JPG, SVG, as well as view/display animations. Some minor editing, such as modifying the lighting or alpha channels, can also be done.

Note
Not implemented for the codimension one case.

Definition at line 2982 of file grid_out.cc.

◆ write_vtk()

template<int dim, int spacedim>
void GridOut::write_vtk ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out 
) const

Write triangulation in VTK format. This function writes a UNSTRUCTURED_GRID file, that contains the following VTK cell types: VTK_HEXAHEDRON, VTK_QUAD, and VTK_LINE, depending on the template dimension.

In three dimensions, this function writes a file that contains

  • VTK_HEXAHEDRON cell types, containing the cell information of the Triangulation
  • VTK_QUAD cell types, containing all boundary faces with non-zero boundary ids, and all faces with non-flat manifold ids
  • VTK_LINE cell types, containing all boundary edges with non-zero boundary ids, and all edges with non-flat manifold ids

In two dimensions:

  • VTK_QUAD cell types, containing the cell information of the Triangulation
  • VTK_LINE cell types, containing all boundary faces with non-zero boundary ids, and all faces with non-flat manifold ids

In one dimension

  • VTK_LINE cell types, containing the cell information of the Triangulation

The output file will contain two CELL_DATA sections, MaterialID and ManifoldID, recording for each VTK cell type the material or boundary id, and the manifold. See the VTK file format documentation for an explanation of the generated output.

The companion GridIn::read_vtk function can be used to read VTK files generated with this method.

Definition at line 3305 of file grid_out.cc.

◆ write_vtu()

template<int dim, int spacedim>
void GridOut::write_vtu ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out 
) const

Write triangulation in VTU format.

Due to the way this function writes data to the output stream, the resulting output files correspond to a faithful representation of the mesh in that all cells are visible for visualization. In general, the data is not in a format that allows reading this file in again through the GridIn class. This is because every vertex of the mesh is duplicated as many times as there are adjacent cells. In other words, every cell has its own, separate set of vertices that are at the same location as the vertices of other cells, but are separately numbered.

In order to create a file that can be read with the GridIn class, the flag GridOutFlags::Vtu::serialize_triangulation must be set to true. In this case, the generated vtu file will contain the triangulation in a xml section which is ignored by general vtu readers.

Definition at line 3544 of file grid_out.cc.

◆ write_mesh_per_processor_as_vtu()

template<int dim, int spacedim>
void GridOut::write_mesh_per_processor_as_vtu ( const Triangulation< dim, spacedim > &  tria,
const std::string &  filename_without_extension,
const bool  view_levels = false,
const bool  include_artificial = false 
) const

Write triangulation in VTU format for each processor, and add a .pvtu file for visualization in VisIt or Paraview that describes the collection of VTU files as all part of the same simulation. The output is in the form filename_without_extension.proc000*.vtu where * is 0,1,...,n_proc-1 and filename_without_extension.pvtu. The input view_levels can be set as true to view each level of a multilevel method. The input include_artificial can be set as true to view the artificial cells for each processor. Each .vtu and .pvtu file will have the attributes subdomain, level_subdomain, level, and proc_writing. The level value can be used to separate the image into the view of the grid on each level of a multilevel method and the proc_writing value can be used to separate the image into each processor's owned and ghost cells. This is accomplished by applying the "warp by scalar" filter in paraview to each of the values. After opening the .pvtu file of a mesh where the input view_levels is set to true, select the "warp by scalar" filter. For the "Scalars" input select proc_writing and for the "Normal" input enter in 1 0 0, then click Apply. Next select the "warp by scalar" filter again. For the "Scalars" input select level and for the "Normal" input enter in 0 1 0, then click Apply. This will give you the following image.

If the view_levels remains at false, thereby only giving the mesh for the active level, it is enough to separate the image into the views written by different processors. This is shown in the following image where the include_artificial input is set as true.

Note: Depending on the size of the mesh you may need to increase the "Scale Factor" input so that each piece does not overlap.

Definition at line 3591 of file grid_out.cc.

◆ write() [1/2]

template<int dim, int spacedim>
void GridOut::write ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out,
const OutputFormat  output_format,
const Mapping< dim, spacedim > *  mapping = nullptr 
) const

Write grid to out according to the given data format. This function simply calls the appropriate write_* function.

Definition at line 5119 of file grid_out.cc.

◆ write() [2/2]

template<int dim, int spacedim>
void GridOut::write ( const Triangulation< dim, spacedim > &  tria,
std::ostream &  out,
const Mapping< dim, spacedim > *  mapping = nullptr 
) const

Write mesh in default format set by ParameterHandler.

Definition at line 5176 of file grid_out.cc.

◆ set_flags() [1/12]

void GridOut::set_flags ( const GridOutFlags::DX flags)

Set flags for DX output

Definition at line 471 of file grid_out.cc.

◆ set_flags() [2/12]

void GridOut::set_flags ( const GridOutFlags::Msh flags)

Set flags for Gmsh output

Definition at line 479 of file grid_out.cc.

◆ set_flags() [3/12]

void GridOut::set_flags ( const GridOutFlags::Ucd flags)

Set flags for UCD output

Definition at line 486 of file grid_out.cc.

◆ set_flags() [4/12]

void GridOut::set_flags ( const GridOutFlags::Gnuplot flags)

Set flags for GNUPLOT output

Definition at line 494 of file grid_out.cc.

◆ set_flags() [5/12]

void GridOut::set_flags ( const GridOutFlags::Eps< 1 > &  flags)

Set flags for EPS output of a one-dimensional triangulation

Definition at line 502 of file grid_out.cc.

◆ set_flags() [6/12]

void GridOut::set_flags ( const GridOutFlags::Eps< 2 > &  flags)

Set flags for EPS output of a two-dimensional triangulation

Definition at line 510 of file grid_out.cc.

◆ set_flags() [7/12]

void GridOut::set_flags ( const GridOutFlags::Eps< 3 > &  flags)

Set flags for EPS output of a three-dimensional triangulation

Definition at line 518 of file grid_out.cc.

◆ set_flags() [8/12]

void GridOut::set_flags ( const GridOutFlags::XFig flags)

Set flags for EPS output of a three-dimensional triangulation

Definition at line 526 of file grid_out.cc.

◆ set_flags() [9/12]

void GridOut::set_flags ( const GridOutFlags::Svg flags)

Set flags for SVG output

Definition at line 533 of file grid_out.cc.

◆ set_flags() [10/12]

void GridOut::set_flags ( const GridOutFlags::MathGL flags)

Set flags for MathGL output

Definition at line 540 of file grid_out.cc.

◆ set_flags() [11/12]

void GridOut::set_flags ( const GridOutFlags::Vtk flags)

Set flags for VTK output

Definition at line 546 of file grid_out.cc.

◆ set_flags() [12/12]

void GridOut::set_flags ( const GridOutFlags::Vtu flags)

Set flags for VTU output

Definition at line 552 of file grid_out.cc.

◆ default_suffix() [1/2]

std::string GridOut::default_suffix ( const OutputFormat  output_format)
static

Provide a function that can tell us which suffix a given output format usually has. For example, it defines the following mappings:

  • OpenDX: .dx
  • gnuplot: .gnuplot
  • ucd: .inp
  • eps: .eps.

Similar mappings are provided for all implemented formats.

Since this function does not need data from this object, it is static and can thus be called without creating an object of this class.

Definition at line 558 of file grid_out.cc.

◆ default_suffix() [2/2]

std::string GridOut::default_suffix ( ) const

Default suffix for the default output format selected through ParameterHandler.

Definition at line 593 of file grid_out.cc.

◆ parse_output_format()

GridOut::OutputFormat GridOut::parse_output_format ( const std::string &  format_name)
static

Return the OutputFormat value corresponding to the given string. If the string does not match any known format, an exception is thrown.

Since this function does not need data from this object, it is static and can thus be called without creating an object of this class. Its main purpose is to allow a program to use any implemented output format without the need to extend the program's parser each time a new format is implemented.

To get a list of presently available format names, e.g. to give it to the ParameterHandler class, use the function get_output_format_names().

Definition at line 601 of file grid_out.cc.

◆ get_output_format_names()

std::string GridOut::get_output_format_names ( )
static

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

Definition at line 644 of file grid_out.cc.

◆ declare_parameters()

void GridOut::declare_parameters ( ParameterHandler param)
static

Declare parameters in ParameterHandler.

Definition at line 651 of file grid_out.cc.

◆ parse_parameters()

void GridOut::parse_parameters ( ParameterHandler param)

Parse parameters of ParameterHandler.

Definition at line 700 of file grid_out.cc.

◆ memory_consumption()

std::size_t GridOut::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 746 of file grid_out.cc.

◆ write_msh_faces() [1/4]

template<int dim, int spacedim>
unsigned int GridOut::write_msh_faces ( const Triangulation< dim, spacedim > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Write the grid information about faces to out. Only those faces are printed which are on the boundary and which have a boundary indicator not equal to zero, since the latter is the default for boundary faces.

Since, in Gmsh, geometric elements are continuously numbered, this function requires a parameter next_element_index providing the next geometric element number. This index should have a numerical value equal to one more than the index previously used to write a geometric element to out.

Returns
The next unused geometric element index.
Warning
next_element_index should be (at least) one larger than the current number of triangulation elements (lines, cells, faces) that have been written to out. Gmsh will not load the saved file correctly if there are repeated indices.

This function unfortunately can not be included in the regular write_msh function, since it needs special treatment for the case dim==1, in which case the face iterators are void*'s and lack the member functions which are called. We would not actually call these functions, but the compiler would complain anyway when compiling the function for dim==1. Bad luck.

Definition at line 3899 of file grid_out.cc.

◆ write_msh_faces() [2/4]

unsigned int GridOut::write_msh_faces ( const Triangulation< 1, 1 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d. Does nothing.

Definition at line 3828 of file grid_out.cc.

◆ write_msh_faces() [3/4]

unsigned int GridOut::write_msh_faces ( const Triangulation< 1, 2 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d, 2sd. Does nothing.

Definition at line 3837 of file grid_out.cc.

◆ write_msh_faces() [4/4]

unsigned int GridOut::write_msh_faces ( const Triangulation< 1, 3 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Definition at line 3845 of file grid_out.cc.

◆ write_msh_lines() [1/6]

template<int dim, int spacedim>
unsigned int GridOut::write_msh_lines ( const Triangulation< dim, spacedim > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Write the grid information about lines to out. Only those lines are printed which are on the boundary and which have a boundary indicator not equal to zero, since the latter is the default for boundary faces.

Since, in Gmsh, geometric elements are continuously numbered, this function requires a parameter next_element_index providing the next geometric element number. This index should have a numerical value equal to one more than the index previously used to write a geometric element to out.

Returns
The next unused geometric element index.
Warning
next_element_index should be (at least) one larger than the current number of triangulation elements (lines, cells, faces) that have been written to out. Gmsh will not load the saved file correctly if there are repeated indices.

This function unfortunately can not be included in the regular write_msh function, since it needs special treatment for the case dim==1 and dim==2, in which case the edge iterators are void*'s and lack the member functions which are called. We would not actually call these functions, but the compiler would complain anyway when compiling the function for dim==1/2. Bad luck.

Definition at line 3938 of file grid_out.cc.

◆ write_msh_lines() [2/6]

unsigned int GridOut::write_msh_lines ( const Triangulation< 1, 1 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d. Does nothing.

Definition at line 3854 of file grid_out.cc.

◆ write_msh_lines() [3/6]

unsigned int GridOut::write_msh_lines ( const Triangulation< 1, 2 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d, 2sd. Does nothing.

Definition at line 3862 of file grid_out.cc.

◆ write_msh_lines() [4/6]

unsigned int GridOut::write_msh_lines ( const Triangulation< 1, 3 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d, 3sd. Does nothing.

Definition at line 3871 of file grid_out.cc.

◆ write_msh_lines() [5/6]

unsigned int GridOut::write_msh_lines ( const Triangulation< 2, 2 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 2d. Does nothing.

Definition at line 3880 of file grid_out.cc.

◆ write_msh_lines() [6/6]

unsigned int GridOut::write_msh_lines ( const Triangulation< 2, 3 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 2d, 3sd. Does nothing.

Definition at line 3888 of file grid_out.cc.

◆ write_ucd_faces() [1/4]

template<int dim, int spacedim>
unsigned int GridOut::write_ucd_faces ( const Triangulation< dim, spacedim > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Write the grid information about faces to out. Only those faces are printed which are on the boundary and which have a boundary indicator not equal to zero, since the latter is the default for boundary faces.

Since (in the UCD format) geometric elements are continuously numbered, this function requires a parameter next_element_index providing the next geometric element number. This index should have a numerical value equal to one more than the index previously used to write a geometric element to out.

Returns
The next unused geometric element index.
Warning
next_element_index should be (at least) one larger than the current number of triangulation elements (lines, cells, faces) that have been written to out. Visualization programs may not load the saved file correctly if there are repeated indices.

This function unfortunately can not be included in the regular write_ucd function, since it needs special treatment for the case dim==1, in which case the face iterators are void*'s and lack the member functions which are called. We would not actually call these functions, but the compiler would complain anyway when compiling the function for dim==1. Bad luck.

Definition at line 4058 of file grid_out.cc.

◆ write_ucd_faces() [2/4]

unsigned int GridOut::write_ucd_faces ( const Triangulation< 1, 1 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d. Does nothing.

Definition at line 3989 of file grid_out.cc.

◆ write_ucd_faces() [3/4]

unsigned int GridOut::write_ucd_faces ( const Triangulation< 1, 2 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d, 2sd. Does nothing.

Definition at line 3997 of file grid_out.cc.

◆ write_ucd_faces() [4/4]

unsigned int GridOut::write_ucd_faces ( const Triangulation< 1, 3 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Definition at line 4005 of file grid_out.cc.

◆ write_ucd_lines() [1/6]

template<int dim, int spacedim>
unsigned int GridOut::write_ucd_lines ( const Triangulation< dim, spacedim > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Write the grid information about lines to out. Only those lines are printed which are on the boundary and which have a boundary indicator not equal to zero, since the latter is the default for boundary lines.

Since (in the UCD format) geometric elements are continuously numbered, this function requires a parameter next_element_index providing the next geometric element number. This index should have a numerical value equal to one more than the index previously used to write a geometric element to out.

Returns
The next unused geometric element index.
Warning
next_element_index should be (at least) one larger than the current number of triangulation elements (lines, cells, faces) that have been written to out. Visualization programs may not load the saved file correctly if there are repeated indices.

This function unfortunately can not be included in the regular write_ucd function, since it needs special treatment for the case dim==1/2, in which case the edge iterators are void*'s and lack the member functions which are called. We would not actually call these functions, but the compiler would complain anyway when compiling the function for dim==1/2. Bad luck.

Definition at line 4100 of file grid_out.cc.

◆ write_ucd_lines() [2/6]

unsigned int GridOut::write_ucd_lines ( const Triangulation< 1, 1 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d. Does nothing.

Definition at line 4013 of file grid_out.cc.

◆ write_ucd_lines() [3/6]

unsigned int GridOut::write_ucd_lines ( const Triangulation< 1, 2 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d, 2sd. Does nothing.

Definition at line 4021 of file grid_out.cc.

◆ write_ucd_lines() [4/6]

unsigned int GridOut::write_ucd_lines ( const Triangulation< 1, 3 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 1d, 3sd. Does nothing.

Definition at line 4030 of file grid_out.cc.

◆ write_ucd_lines() [5/6]

unsigned int GridOut::write_ucd_lines ( const Triangulation< 2, 2 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 2d. Does nothing.

Definition at line 4039 of file grid_out.cc.

◆ write_ucd_lines() [6/6]

unsigned int GridOut::write_ucd_lines ( const Triangulation< 2, 3 > &  tria,
const unsigned int  next_element_index,
std::ostream &  out 
) const
private

Declaration of the specialization of above function for 2d, 3sd. Does nothing.

Definition at line 4047 of file grid_out.cc.

◆ n_boundary_faces() [1/4]

template<int dim, int spacedim>
unsigned int GridOut::n_boundary_faces ( const Triangulation< dim, spacedim > &  tria) const
private

Return the number of faces in the triangulation which have a boundary indicator not equal to zero. Only these faces are explicitly printed in the write_* functions; all faces with indicator numbers::internal_face_boundary_id are interior ones and an indicator with value zero for faces at the boundary are considered default.

This function always returns an empty list in one dimension.

The reason for this function is the same as for write_ucd_faces(). See there for more information.

Definition at line 3779 of file grid_out.cc.

◆ n_boundary_faces() [2/4]

unsigned int GridOut::n_boundary_faces ( const Triangulation< 1, 1 > &  tria) const
private

Declaration of the specialization of above function for 1d. Simply returns zero.

Definition at line 3727 of file grid_out.cc.

◆ n_boundary_faces() [3/4]

unsigned int GridOut::n_boundary_faces ( const Triangulation< 1, 2 > &  tria) const
private

Declaration of the specialization of above function for 1d, 2sd. Simply returns zero.

Definition at line 3740 of file grid_out.cc.

◆ n_boundary_faces() [4/4]

unsigned int GridOut::n_boundary_faces ( const Triangulation< 1, 3 > &  tria) const
private

Declaration of the specialization of above function for 1d, 3sd. Simply returns zero.

Definition at line 3752 of file grid_out.cc.

◆ n_boundary_lines() [1/6]

template<int dim, int spacedim>
unsigned int GridOut::n_boundary_lines ( const Triangulation< dim, spacedim > &  tria) const
private

Return the number of lines in the triangulation which have a boundary indicator not equal to zero. Only these lines are explicitly printed in the write_* functions; all lines with indicator numbers::internal_face_boundary_id are interior ones and an indicator with value zero for faces at the boundary are considered default.

This function always returns an empty list in one and two dimensions.

The reason for this function is the same as for write_ucd_faces(). See there for more information.

Definition at line 3795 of file grid_out.cc.

◆ n_boundary_lines() [2/6]

unsigned int GridOut::n_boundary_lines ( const Triangulation< 1, 1 > &  tria) const
private

Declaration of the specialization of above function for 1d. Simply returns zero.

Definition at line 3733 of file grid_out.cc.

◆ n_boundary_lines() [3/6]

unsigned int GridOut::n_boundary_lines ( const Triangulation< 1, 2 > &  tria) const
private

Declaration of the specialization of above function for 1d, 2sd. Simply returns zero.

Definition at line 3746 of file grid_out.cc.

◆ n_boundary_lines() [4/6]

unsigned int GridOut::n_boundary_lines ( const Triangulation< 1, 3 > &  tria) const
private

Declaration of the specialization of above function for 1d, 3sd. Simply returns zero.

Definition at line 3758 of file grid_out.cc.

◆ n_boundary_lines() [5/6]

unsigned int GridOut::n_boundary_lines ( const Triangulation< 2, 2 > &  tria) const
private

Declaration of the specialization of above function for 2d. Simply returns zero.

Definition at line 3764 of file grid_out.cc.

◆ n_boundary_lines() [6/6]

unsigned int GridOut::n_boundary_lines ( const Triangulation< 2, 3 > &  tria) const
private

Declaration of the specialization of above function for 2d, 3sd. Simply returns zero.

Definition at line 3770 of file grid_out.cc.

◆ write_dx() [2/4]

template<>
void GridOut::write_dx ( const Triangulation< 1 > &  ,
std::ostream &   
) const

Definition at line 758 of file grid_out.cc.

◆ write_dx() [3/4]

template<>
void GridOut::write_dx ( const Triangulation< 1, 2 > &  ,
std::ostream &   
) const

Definition at line 765 of file grid_out.cc.

◆ write_dx() [4/4]

template<>
void GridOut::write_dx ( const Triangulation< 1, 3 > &  ,
std::ostream &   
) const

Definition at line 772 of file grid_out.cc.

◆ write_xfig() [2/2]

template<>
void GridOut::write_xfig ( const Triangulation< 2 > &  tria,
std::ostream &  out,
const Mapping< 2 > *   
) const

Definition at line 1254 of file grid_out.cc.

◆ write_mathgl() [2/2]

template<>
void GridOut::write_mathgl ( const Triangulation< 1 > &  ,
std::ostream &   
) const

Definition at line 2972 of file grid_out.cc.

Member Data Documentation

◆ default_format

OutputFormat GridOut::default_format
private

The default output format, set by a ParameterHandler.

Definition at line 1549 of file grid_out.h.

◆ dx_flags

GridOutFlags::DX GridOut::dx_flags
private

Flags for OpenDX output.

Definition at line 1554 of file grid_out.h.

◆ msh_flags

GridOutFlags::Msh GridOut::msh_flags
private

Flags for Gmsh output. Can be changed by using the set_flags(const GridOutFlags::Msh&) function.

Definition at line 1560 of file grid_out.h.

◆ ucd_flags

GridOutFlags::Ucd GridOut::ucd_flags
private

Flags for UCD output. Can be changed by using the set_flags(const GridOutFlags::Ucd&) function.

Definition at line 1566 of file grid_out.h.

◆ gnuplot_flags

GridOutFlags::Gnuplot GridOut::gnuplot_flags
private

Flags to be used upon output of GNUPLOT data. Can be changed by using the set_flags(const GridOutFlags::Gnuplot&) function.

Definition at line 1572 of file grid_out.h.

◆ eps_flags_1

GridOutFlags::Eps<1> GridOut::eps_flags_1
private

Flags to be used upon output of EPS data in one space dimension. Can be changed by using the set_flags(const GridOutFlags::Eps<1>&) function.

Definition at line 1578 of file grid_out.h.

◆ eps_flags_2

GridOutFlags::Eps<2> GridOut::eps_flags_2
private

Flags to be used upon output of EPS data in two space dimensions. Can be changed by using the set_flags function.

Definition at line 1584 of file grid_out.h.

◆ eps_flags_3

GridOutFlags::Eps<3> GridOut::eps_flags_3
private

Flags to be used upon output of EPS data in three space dimensions. Can be changed by using the set_flags function.

Definition at line 1590 of file grid_out.h.

◆ xfig_flags

GridOutFlags::XFig GridOut::xfig_flags
private

Flags used for XFig output.

Definition at line 1595 of file grid_out.h.

◆ svg_flags

GridOutFlags::Svg GridOut::svg_flags
private

Flags used for Svg output.

Definition at line 1600 of file grid_out.h.

◆ mathgl_flags

GridOutFlags::MathGL GridOut::mathgl_flags
private

Flags for MathGL output.

Definition at line 1605 of file grid_out.h.

◆ vtk_flags

GridOutFlags::Vtk GridOut::vtk_flags
private

Flags for VTK output.

Definition at line 1610 of file grid_out.h.

◆ vtu_flags

GridOutFlags::Vtu GridOut::vtu_flags
private

Flags for VTU output.

Definition at line 1615 of file grid_out.h.


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