Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d0/dfe/namespacedolfinx_1_1io_1_1xdmf__mesh.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
Functions
dolfinx::io::xdmf_mesh Namespace Reference

Low-level methods for reading XDMF files. More...

Functions

template<std::floating_point U>
void add_mesh (MPI_Comm comm, pugi::xml_node &xml_node, hid_t h5_id, const mesh::Mesh< U > &mesh, const std::string &path_prefix)
 Add Mesh to xml node.
 
template<std::floating_point U>
void add_topology_data (MPI_Comm comm, pugi::xml_node &xml_node, hid_t h5_id, std::string path_prefix, const mesh::Topology &topology, const mesh::Geometry< U > &geometry, int cell_dim, std::span< const std::int32_t > entities)
 Add Topology xml node.
 
template<std::floating_point U>
void add_geometry_data (MPI_Comm comm, pugi::xml_node &xml_node, hid_t h5_id, std::string path_prefix, const mesh::Geometry< U > &geometry)
 Add Geometry xml node.
 
std::pair< std::variant< std::vector< float >, std::vector< double > >, std::array< std::size_t, 2 > > read_geometry_data (MPI_Comm comm, hid_t h5_id, const pugi::xml_node &node)
 Read geometry (coordinate) data.
 
std::pair< std::vector< std::int64_t >, std::array< std::size_t, 2 > > read_topology_data (MPI_Comm comm, hid_t h5_id, const pugi::xml_node &node)
 Read topology (cell connectivity) data.
 
template<typename T , std::floating_point U>
void add_meshtags (MPI_Comm comm, const mesh::MeshTags< T > &meshtags, const mesh::Geometry< U > &geometry, pugi::xml_node &xml_node, hid_t h5_id, const std::string &name)
 Add mesh tags to XDMF file.
 

Detailed Description

Low-level methods for reading XDMF files.

Function Documentation

◆ add_mesh()

template<std::floating_point U>
void add_mesh ( MPI_Comm  comm,
pugi::xml_node &  xml_node,
hid_t  h5_id,
const mesh::Mesh< U > &  mesh,
const std::string &  path_prefix 
)

Add Mesh to xml node.

Creates new Grid with Topology and Geometry xml nodes for mesh. In HDF file data is stored under path prefix.

◆ add_topology_data()

template<std::floating_point U>
void add_topology_data ( MPI_Comm  comm,
pugi::xml_node &  xml_node,
hid_t  h5_id,
std::string  path_prefix,
const mesh::Topology topology,
const mesh::Geometry< U > &  geometry,
int  cell_dim,
std::span< const std::int32_t >  entities 
)

Add Topology xml node.

Parameters
[in]comm
[in]xml_node
[in]h5_id
[in]path_prefix
[in]topology
[in]geometry
[in]cell_dimDimension of mesh entities to save
[in]entitiesLocal-to-process indices of mesh entities whose topology will be saved. This is used to save subsets of Mesh.

◆ read_geometry_data()

std::pair< std::variant< std::vector< float >, std::vector< double > >, std::array< std::size_t, 2 > > read_geometry_data ( MPI_Comm  comm,
hid_t  h5_id,
const pugi::xml_node &  node 
)

Read geometry (coordinate) data.

Returns
The coordinates of each 'node'. The returned data is (0) an array holding the coordinates (row-major storage) and (1) the shape of the coordinate array. The shape is (num_nodes, geometric dimension).

◆ read_topology_data()

std::pair< std::vector< std::int64_t >, std::array< std::size_t, 2 > > read_topology_data ( MPI_Comm  comm,
hid_t  h5_id,
const pugi::xml_node &  node 
)

Read topology (cell connectivity) data.

Returns
Mesh topology in DOLFINx ordering, where data row i lists the 'nodes' of cell i. The returned data is (0) an array holding the topology data (row-major storage) and (1) the shape of the topology array. The shape is (num_cells, num_nodes_per_cell)