Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d9/d9e/classdolfinx_1_1io_1_1XDMFFile.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Attributes | List of all members
XDMFFile Class Reference

Read and write mesh::Mesh, fem::Function and other objects in XDMF. More...

#include <XDMFFile.h>

Public Types

enum class  Encoding { HDF5 , ASCII }
 File encoding type.
 

Public Member Functions

 XDMFFile (MPI_Comm comm, const std::filesystem::path &filename, std::string file_mode, Encoding encoding=default_encoding)
 Constructor.
 
 XDMFFile (XDMFFile &&)=default
 Move constructor.
 
 ~XDMFFile ()
 Destructor.
 
void close ()
 Close the file.
 
template<std::floating_point U>
void write_mesh (const mesh::Mesh< U > &mesh, std::string xpath="/Xdmf/Domain")
 Save Mesh.
 
void write_geometry (const mesh::Geometry< double > &geometry, std::string name, std::string xpath="/Xdmf/Domain")
 Save Geometry.
 
mesh::Mesh< double > read_mesh (const fem::CoordinateElement< double > &element, mesh::GhostMode mode, std::string name, std::string xpath="/Xdmf/Domain") const
 Read in Mesh.
 
std::pair< std::vector< std::int64_t >, std::array< std::size_t, 2 > > read_topology_data (std::string name, std::string xpath="/Xdmf/Domain") const
 Read Topology data for Mesh.
 
std::pair< std::variant< std::vector< float >, std::vector< double > >, std::array< std::size_t, 2 > > read_geometry_data (std::string name, std::string xpath="/Xdmf/Domain") const
 Read Geometry data for Mesh.
 
std::pair< mesh::CellType, int > read_cell_type (std::string grid_name, std::string xpath="/Xdmf/Domain")
 Read information about cell type.
 
template<dolfinx::scalar T, std::floating_point U = scalar_value_type_t<T>>
void write_function (const fem::Function< T, U > &u, double t, std::string mesh_xpath="/Xdmf/Domain/Grid[@GridType='Uniform'][1]")
 Write a fem::Function to file.
 
template<std::floating_point T>
void write_meshtags (const mesh::MeshTags< std::int32_t > &meshtags, const mesh::Geometry< T > &x, std::string geometry_xpath, std::string xpath="/Xdmf/Domain")
 Write MeshTags.
 
mesh::MeshTags< std::int32_t > read_meshtags (const mesh::Mesh< double > &mesh, std::string name, std::string xpath="/Xdmf/Domain")
 Read MeshTags.
 
void write_information (std::string name, std::string value, std::string xpath="/Xdmf/Domain/")
 Write Information.
 
std::string read_information (std::string name, std::string xpath="/Xdmf/Domain/")
 Read Information.
 
MPI_Comm comm () const
 Get the MPI communicator.
 

Static Public Attributes

static const Encoding default_encoding = Encoding::HDF5
 Default encoding type.
 

Detailed Description

Read and write mesh::Mesh, fem::Function and other objects in XDMF.

This class supports the output of meshes and functions in XDMF (http://www.xdmf.org) format. It creates an XML file that describes the data and points to a HDF5 file that stores the actual problem data. Output of data in parallel is supported.

XDMF is not suitable for higher order geometries, as their currently only supports 1st and 2nd order geometries.

Member Function Documentation

◆ close()

void close ( )

Close the file.

This closes open underlying HDF5 file. In ASCII mode the XML file is closed each time it is written to or read from, so close() has no effect.

◆ comm()

MPI_Comm comm ( ) const

Get the MPI communicator.

Returns
The MPI communicator for the file object

◆ read_cell_type()

std::pair< mesh::CellType, int > read_cell_type ( std::string  grid_name,
std::string  xpath = "/Xdmf/Domain" 
)

Read information about cell type.

Parameters
[in]grid_nameName of Grid for which cell type is needed
[in]xpathXPath where Grid is stored

◆ read_geometry_data()

std::pair< std::variant< std::vector< float >, std::vector< double > >, std::array< std::size_t, 2 > > read_geometry_data ( std::string  name,
std::string  xpath = "/Xdmf/Domain" 
) const

Read Geometry data for Mesh.

Parameters
[in]nameName of the mesh (Grid)
[in]xpathXPath where Mesh Grid data is located
Returns
points on each process

◆ read_information()

std::string read_information ( std::string  name,
std::string  xpath = "/Xdmf/Domain/" 
)

Read Information.

Parameters
[in]name
[in]xpathXPath where Information is stored in file

◆ read_mesh()

mesh::Mesh< double > read_mesh ( const fem::CoordinateElement< double > &  element,
mesh::GhostMode  mode,
std::string  name,
std::string  xpath = "/Xdmf/Domain" 
) const

Read in Mesh.

Parameters
[in]elementElement that describes the geometry of a cell
[in]modeThe type of ghosting/halo to use for the mesh when distributed in parallel
[in]name
[in]xpathXPath where Mesh Grid is located
Returns
A Mesh distributed on the same communicator as the XDMFFile

◆ read_meshtags()

mesh::MeshTags< std::int32_t > read_meshtags ( const mesh::Mesh< double > &  mesh,
std::string  name,
std::string  xpath = "/Xdmf/Domain" 
)

Read MeshTags.

Parameters
[in]meshThe Mesh that the data is defined on
[in]name
[in]xpathXPath where MeshTags Grid is stored in file

◆ read_topology_data()

std::pair< std::vector< std::int64_t >, std::array< std::size_t, 2 > > read_topology_data ( std::string  name,
std::string  xpath = "/Xdmf/Domain" 
) const

Read Topology data for Mesh.

Parameters
[in]nameName of the mesh (Grid)
[in]xpathXPath where Mesh Grid data is located
Returns
(Cell type, degree), and cells topology (global node indexing)

◆ write_function()

template<dolfinx::scalar T, std::floating_point U>
void write_function ( const fem::Function< T, U > &  u,
double  t,
std::string  mesh_xpath = "/Xdmf/Domain/Grid[@GridType='Uniform'][1]" 
)

Write a fem::Function to file.

Precondition
The fem::Function u must be (i) a lowest-order (P0) discontinuous Lagrange element or (ii) a continuous Lagrange element where the element 'nodes' are the same as the nodes of its mesh::Mesh. Otherwise an exception is raised.
Note
User interpolation to a suitable Lagrange space may be required to satisfy the precondition on u. The VTX output (io::VTXWriter) format is recommended over XDMF for discontinuous and/or high-order spaces.
Parameters
[in]uFunction to write to file.
[in]tTime stamp to associate with u.
[in]mesh_xpathXPath for a Grid under which u will be inserted.

◆ write_geometry()

void write_geometry ( const mesh::Geometry< double > &  geometry,
std::string  name,
std::string  xpath = "/Xdmf/Domain" 
)

Save Geometry.

Parameters
[in]geometry
[in]name
[in]xpathXPath of a node where Geometry will be inserted

◆ write_information()

void write_information ( std::string  name,
std::string  value,
std::string  xpath = "/Xdmf/Domain/" 
)

Write Information.

Parameters
[in]name
[in]valueString to store into Information tag
[in]xpathXPath where Information will be inserted

◆ write_mesh()

template<std::floating_point U>
void write_mesh ( const mesh::Mesh< U > &  mesh,
std::string  xpath = "/Xdmf/Domain" 
)

Save Mesh.

Parameters
[in]mesh
[in]xpathXPath where Mesh Grid will be written

◆ write_meshtags()

template<std::floating_point T>
void write_meshtags ( const mesh::MeshTags< std::int32_t > &  meshtags,
const mesh::Geometry< T > &  x,
std::string  geometry_xpath,
std::string  xpath = "/Xdmf/Domain" 
)

Write MeshTags.

Parameters
[in]meshtags
[in]xMesh geometry
[in]geometry_xpathXPath where Geometry is already stored in file
[in]xpathXPath where MeshTags Grid will be inserted

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