dolfinx.io
Tools for file input/output (IO).
Functions
|
Given a set of mesh entities and values, distribute them to the process that owns the entity. |
Classes
|
Interface to VTK files. |
|
|
|
Writer for VTX files, using ADIOS2 to create the files. |
|
- class dolfinx.io.VTKFile(self, comm: MPICommWrapper, filename: str | os.PathLike, mode: str)[source]
Bases:
VTKFile
Interface to VTK files.
VTK supports arbitrary order Lagrange finite elements for the geometry description. XDMF is the preferred format for geometry order <= 2.
- class dolfinx.io.VTXMeshPolicy(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)
Bases:
Enum
- reuse = 1
- update = 0
- class dolfinx.io.VTXWriter(comm: Comm, filename: str | Path, output: Mesh | Function | list[Function] | tuple[Function], engine: str = 'BPFile', mesh_policy: VTXMeshPolicy = VTXMeshPolicy.update)[source]
Bases:
object
Writer for VTX files, using ADIOS2 to create the files.
VTX supports arbitrary order Lagrange finite elements for the geometry description and arbitrary order (discontinuous) Lagrange finite elements for Functions.
The files can be viewed using Paraview.
Initialize a writer for outputting data in the VTX format.
- Parameters:
comm – The MPI communicator
filename – The output filename
output – The data to output. Either a mesh, a single (discontinuous) Lagrange Function or list of (discontinuous) Lagrange Functions.
engine – ADIOS2 engine to use for output. See ADIOS2 documentation for options.
mesh_policy – Controls if the mesh is written to file at the first time step only when a
Function
is written to file, or is re-written (updated) at each time step. Has an effect only forFunction
output.
Note
All Functions for output must share the same mesh and have the same element type.
- class dolfinx.io.XDMFFile(self, comm: MPICommWrapper, filename: str | os.PathLike, file_mode: str, encoding: dolfinx.cpp.io.XDMFFile.Encoding = Encoding.HDF5)[source]
Bases:
XDMFFile
- read_mesh(ghost_mode=GhostMode.shared_facet, name='mesh', xpath='/Xdmf/Domain') Mesh [source]
Read mesh data from file.
- read_meshtags(self, mesh: dolfinx.cpp.mesh.Mesh_float64, name: str, xpath: str = '/Xdmf/Domain') dolfinx.cpp.mesh.MeshTags_int32 [source]
- write_function(u: Function, t: float = 0.0, mesh_xpath="/Xdmf/Domain/Grid[@GridType='Uniform'][1]")[source]
Write function to file for a given time.
Note
Function is interpolated onto the mesh nodes, as a Nth order Lagrange function, where N is the order of the coordinate map. If the Function is a cell-wise constant, it is saved as a cell-wise constant.
- Parameters:
u – Function to write to file.
t – Time associated with Function output.
mesh_xpath – Path to mesh associated with the Function in the XDMFFile.
- dolfinx.io.distribute_entity_data(mesh: Mesh, entity_dim: int, entities: ndarray[Any, dtype[int64]], values: ndarray) tuple[ndarray[Any, dtype[int64]], ndarray] [source]
Given a set of mesh entities and values, distribute them to the process that owns the entity.
The entities are described by the global vertex indices of the mesh. These entity indices are using the original input ordering.
- Returns:
Entities owned by the process (and their local entity-to-vertex indices) and the corresponding values.