dolfinx.io

Tools for file input/output (IO).

class dolfinx.io.FidesMeshPolicy(self: dolfinx.cpp.io.FidesMeshPolicy, value: int)

Bases: pybind11_object

Members:

update

reuse

property name
reuse = <FidesMeshPolicy.reuse: 1>
update = <FidesMeshPolicy.update: 0>
property value
class dolfinx.io.FidesWriter(comm: ~mpi4py.MPI.Comm, filename: ~typing.Union[str, ~pathlib.Path], output: ~typing.Union[~dolfinx.mesh.Mesh, ~typing.List[~dolfinx.fem.function.Function], ~dolfinx.fem.function.Function], engine: ~typing.Optional[str] = 'BPFile', mesh_policy: ~typing.Optional[~dolfinx.cpp.io.FidesMeshPolicy] = <FidesMeshPolicy.update: 0>)[source]

Bases: object

Writer for Fides files, using ADIOS2 to create the files.

Fides (https://fides.readthedocs.io/) supports first order Lagrange finite elements for the geometry description and first order Lagrange finite elements for functions. All functions have to be of the same element family and same order.

The files can be displayed by Paraview.

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 – MPI communicator.

  • filename – Output filename.

  • output – Data to output. Either a mesh, a single degree one Lagrange function or list of degree one 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 for Function output.

close()[source]
write(t: float)[source]
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.

write_function(u: Union[List[Function], Function], t: float = 0.0) None[source]

Write a single function or a list of functions to file for a given time (default 0.0)

write_mesh(mesh: Mesh, t: float = 0.0) None[source]

Write mesh to file for a given time (default 0.0)

class dolfinx.io.VTXWriter(comm: Comm, filename: Union[str, Path], output: Union[Mesh, Function, List[Function]], engine: Optional[str] = 'BPFile')[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 displayed by 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.

Note

All Functions for output must share the same mesh and have the same element type.

close()[source]
write(t: float)[source]
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_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.

write_mesh(mesh: Mesh, xpath: str = '/Xdmf/Domain') None[source]

Write mesh to file

write_meshtags(tags: MeshTags, x: Union[Geometry_float32, Geometry_float64], geometry_xpath: str = '/Xdmf/Domain/Grid/Geometry', xpath: str = '/Xdmf/Domain') None[source]

Write mesh tags to file

dolfinx.io.distribute_entity_data(mesh: Mesh, entity_dim: int, entities: ndarray[Any, dtype[int64]], values: ndarray[Any, dtype[int32]]) Tuple[ndarray[Any, dtype[int64]], ndarray[Any, dtype[int32]]][source]