dolfin.io

IO module for input data, post-processing and checkpointing

Classes

HDF5File(mpi_comm, filename, mode)

Interface to HDF5 files

VTKFile(filename)

Interface to VTK files VTK supports arbitrary order Lagrangian finite elements for the geometry description.

XDMFFile(mpi_comm, filename[, encoding])

Interface to XDMF files This format is preferred on lower order geometries and for DG and RT function spaces.

class dolfin.io.HDF5File(mpi_comm, filename: str, mode: str)[source]

Bases: object

Interface to HDF5 files

Open HDF5 file

Parameters
  • mpi_comm – The MPI communicator

  • filename – Name of the file

  • mode – File opening mode, which can be write (w), read (r) or append (a)

close() → None[source]

Close file

get_mpi_atomicity() → bool[source]

Get atomicity of the HDF5 file

read_function(V, name: str)[source]

Read finite element Function from file

Parameters
  • V – Function space of saved function.

  • name – Name of function as saved into HDF file.

Returns

Function read from file

Return type

dolfin.function.function.Function

Note

Parameter V: Function space must be the same as saved function space except for ordering of mesh entities.

read_mf_double(mesh, name: str = '')[source]

Read MeshFunction of type float

read_mvc_size_t(mesh, name: str = '')[source]

Read MeshValueCollection of type size_t

read_vector(mpi_comm, data_path: str, use_partition_from_file: bool)[source]

Read Vector

set_mpi_atomicity(atomicity: bool) → None[source]

Set atomicity of the HDF5 file

write(o, name, t=None) → None[source]

Write object to file

class dolfin.io.XDMFFile(mpi_comm, filename: str, encoding=Encoding.HDF5)[source]

Bases: object

Interface to XDMF files This format is preferred on lower order geometries and for DG and RT function spaces. XDMF also allows for checkpointing of solutions and has parallel support.

Open XDMF file

Parameters
  • mpi_comm – The MPI communicator

  • filename – Name of the file

  • encoding – Encoding used for ‘heavy’ data when writing/appending

class Encoding(self: dolfin.cpp.io.XDMFFile.Encoding, arg0: int) → None

Bases: pybind11_builtins.pybind11_object

Members:

HDF5

ASCII

property name

handle) -> str

Type

(self

close() → None[source]

Close file

read_checkpoint(V, name: str, counter: int = -1) → dolfin.function.Function[source]

Read finite element Function from checkpointing format

Parameters
  • V – Function space of saved function.

  • name – Name of function as saved into XDMF file.

  • counter (optional) – Position of function in the file within functions of the same name. Counter is used to read function saved as time series. To get last saved function use counter=-1, or counter=-2 for one before last, etc.

Note

Parameter V: Function space must be the same as saved function space except for ordering of mesh entities.

Returns

The finite element Function read from checkpoint file

Return type

dolfin.function.function.Function

read_mesh_data(mpi_comm) → Tuple[dolfin.cpp.mesh.CellType, numpy.ndarray, numpy.ndarray, List[int]][source]

Read in mesh data

Parameters

mpi_comm – MPI communicator

Returns

  • cell_type – Cell type

  • points – Geometric points on each process

  • cells – Topological cells with global vertex indexing

  • global_cell_indices – List of global cell indices

read_mf_double(mesh, name: str = '')[source]

Read MeshFunction of type double

read_mf_int(mesh, name: str = '')[source]

Read MeshFunction of type int

read_mf_size_t(mesh, name: str = '')[source]

Read MeshFunction of type size_t

read_mvc_bool(mesh, name: str = '')[source]

Read MeshValueCollection of type bool

read_mvc_double(mesh, name: str = '')[source]

Read MeshValueCollection of type float

read_mvc_int(mesh, name: str = '')[source]

Read MeshValueCollection of type int

read_mvc_size_t(mesh, name: str = '')[source]

Read MeshValueCollection of type size_t

write(o, t=None) → None[source]

Write object to file

Parameters
  • o – The object to write to file

  • t – The time stamp

write_checkpoint(u, name: str, time_step: float = 0.0) → None[source]

Write finite element Function in checkpointing format