dolfinx.io
Tools for input/output (IO).
- class dolfinx.io.FidesWriter(comm: Comm, filename: str, output: Union[Mesh, List[Function], Function])[source]
Bases:
FidesWriter
Interface to Fides file format.
Fides supports first order Lagrange finite elements for the geometry descriptionand first order Lagrange finite elements for functions. All functions has to be of the same element family and same order.
The files can be displayed by Paraview. The storage backend uses ADIOS2.
Initialize a writer for outputting a mesh, a single Lagrange function or list of Lagrange functions sharing the same element family and degree
- Parameters
comm – The MPI communicator
filename – The output filename
output – The data to output. Either a mesh, a single first order Lagrange function or list of first order Lagrange functions.
- class dolfinx.io.VTKFile(self: dolfinx.cpp.io.VTKFile, comm: MPICommWrapper, filename: 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.VTXWriter(comm: Comm, filename: str, output: Union[Mesh, List[Function], Function])[source]
Bases:
VTXWriter
Interface to VTK files for ADIOS2
VTX supports arbitrary order Lagrange finite elements for the geometry description and arbitrary order (discontinuous) Lagrange finite elements for Functions.
The files can be displayed by Paraview. The storage backend uses ADIOS2.
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.
Note
All Functions for output must share the same mesh and have the same element type.
- class dolfinx.io.XDMFFile(self: dolfinx.cpp.io.XDMFFile, comm: MPICommWrapper, filename: os.PathLike, file_mode: str, encoding: dolfinx.cpp.io.XDMFFile.Encoding = <Encoding.HDF5: 0>)[source]
Bases:
XDMFFile
- read_mesh(ghost_mode=<GhostMode.shared_facet: 1>, name='mesh', xpath='/Xdmf/Domain') Mesh [source]
Read mesh data from file
- read_meshtags(self: dolfinx.cpp.io.XDMFFile, mesh: dolfinx.cpp.mesh.Mesh, name: str, xpath: str = '/Xdmf/Domain') dolfinx.cpp.mesh.MeshTags_int32 [source]
- write_function(*args, **kwargs)[source]
Overloaded function.
write_function(self: dolfinx.cpp.io.XDMFFile, function: dolfinx.cpp.fem.Function_float64, t: float, mesh_xpath: str) -> None
write_function(self: dolfinx.cpp.io.XDMFFile, function: dolfinx.cpp.fem.Function_complex128, t: float, mesh_xpath: str) -> None
- dolfinx.io.distribute_entity_data(mesh: dolfinx.cpp.mesh.Mesh, entity_dim: int, entities: numpy.ndarray[numpy.int64], values: numpy.ndarray[numpy.int32]) Tuple[numpy.ndarray[numpy.int32], numpy.ndarray[numpy.int32]]