Read and write mesh::Mesh, fem::Function and other objects in XDMF.
More...
#include <XDMFFile.h>
|
enum class | Encoding : std::int8_t { HDF5
, ASCII
} |
| File encoding type.
|
|
| XDMFFile (MPI_Comm comm, const std::filesystem::path &filename, const std::string &file_mode, Encoding encoding=Encoding::HDF5) |
| Constructor.
|
| XDMFFile (XDMFFile &&)=default |
| Move constructor.
|
| ~XDMFFile () |
| Destructor.
|
void | close () |
template<std::floating_point U> |
void | write_mesh (const mesh::Mesh< U > &mesh, const std::string &xpath="/Xdmf/Domain") |
void | write_geometry (const mesh::Geometry< double > &geometry, const std::string &name, const std::string &xpath="/Xdmf/Domain") |
mesh::Mesh< double > | read_mesh (const fem::CoordinateElement< double > &element, mesh::GhostMode mode, const std::string &name, const std::string &xpath="/Xdmf/Domain") const |
std::pair< std::vector< std::int64_t >, std::array< std::size_t, 2 > > | read_topology_data (const std::string &name, const std::string &xpath="/Xdmf/Domain") const |
std::pair< std::variant< std::vector< float >, std::vector< double > >, std::array< std::size_t, 2 > > | read_geometry_data (const std::string &name, const std::string &xpath="/Xdmf/Domain") const |
std::pair< mesh::CellType, int > | read_cell_type (const std::string &grid_name, const std::string &xpath="/Xdmf/Domain") |
template<dolfinx::scalar T, std::floating_point U = scalar_value_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, const std::string &geometry_xpath, const std::string &xpath="/Xdmf/Domain") |
mesh::MeshTags< std::int32_t > | read_meshtags (const mesh::Mesh< double > &mesh, const std::string &name, std::optional< std::string > attribute_name, const std::string &xpath="/Xdmf/Domain") |
void | write_information (const std::string &name, const std::string &value, const std::string &xpath="/Xdmf/Domain/") |
std::string | read_information (const std::string &name, const std::string &xpath="/Xdmf/Domain/") |
MPI_Comm | comm () const |
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; it only supports 1st- and 2nd-order geometries.
◆ 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()
Get the MPI communicator
- Returns
- The MPI communicator for the file object
◆ read_cell_type()
std::pair< mesh::CellType, int > read_cell_type |
( |
const std::string & | grid_name, |
|
|
const std::string & | xpath = "/Xdmf/Domain" ) |
Read information about cell type
- Parameters
-
[in] | grid_name | Name of Grid for which cell type is needed |
[in] | xpath | XPath 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 |
( |
const std::string & | name, |
|
|
const std::string & | xpath = "/Xdmf/Domain" ) const |
Read Geometry data for Mesh
- Parameters
-
[in] | name | Name of the mesh (Grid) |
[in] | xpath | XPath where Mesh Grid data is located |
- Returns
- points on each process
◆ read_information()
std::string read_information |
( |
const std::string & | name, |
|
|
const std::string & | xpath = "/Xdmf/Domain/" ) |
Read Information
- Parameters
-
[in] | name | |
[in] | xpath | XPath where Information is stored in file |
◆ read_mesh()
Read Mesh
- Parameters
-
[in] | element | Element that describes the geometry of a cell |
[in] | mode | The type of ghosting/halo to use for the mesh when distributed in parallel |
[in] | name | |
[in] | xpath | XPath 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, |
|
|
const std::string & | name, |
|
|
std::optional< std::string > | attribute_name, |
|
|
const std::string & | xpath = "/Xdmf/Domain" ) |
Read MeshTags
- Parameters
-
[in] | mesh | Mesh that the input data is defined on |
[in] | name | Name of the grid node in the xml-scheme of the XDMF-file. E.g. "Material" in Grid Name="Material" GridType="Uniform" |
[in] | attribute_name | Name of the attribute to read |
[in] | xpath | XPath 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 |
( |
const std::string & | name, |
|
|
const std::string & | xpath = "/Xdmf/Domain" ) const |
Read Topology data for Mesh
- Parameters
-
[in] | name | Name of the mesh (Grid) |
[in] | xpath | XPath where Mesh Grid data is located |
- Returns
- (Cell type, degree), and cells topology (global node indexing)
◆ write_function()
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] | u | Function to write to file. |
[in] | t | Time stamp to associate with u. |
[in] | mesh_xpath | XPath for a Grid under which u will be inserted. |
◆ write_geometry()
void write_geometry |
( |
const mesh::Geometry< double > & | geometry, |
|
|
const std::string & | name, |
|
|
const std::string & | xpath = "/Xdmf/Domain" ) |
Save Geometry
- Parameters
-
[in] | geometry | |
[in] | name | |
[in] | xpath | XPath of a node where Geometry will be inserted |
◆ write_information()
void write_information |
( |
const std::string & | name, |
|
|
const std::string & | value, |
|
|
const std::string & | xpath = "/Xdmf/Domain/" ) |
Write Information
- Parameters
-
[in] | name | |
[in] | value | String to store into Information tag |
[in] | xpath | XPath where Information will be inserted |
◆ write_mesh()
template<std::floating_point U>
void write_mesh |
( |
const mesh::Mesh< U > & | mesh, |
|
|
const std::string & | xpath = "/Xdmf/Domain" ) |
Save Mesh
- Parameters
-
[in] | mesh | |
[in] | xpath | XPath 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, |
|
|
const std::string & | geometry_xpath, |
|
|
const std::string & | xpath = "/Xdmf/Domain" ) |
Write MeshTags
- Parameters
-
[in] | meshtags | |
[in] | x | Mesh geometry |
[in] | geometry_xpath | XPath where Geometry is already stored in file |
[in] | xpath | XPath where MeshTags Grid will be inserted |
The documentation for this class was generated from the following files:
- /__w/dolfinx/dolfinx/cpp/dolfinx/io/XDMFFile.h
- /__w/dolfinx/dolfinx/cpp/dolfinx/io/XDMFFile.cpp