XDMFFile

class dolfin.cpp.io.XDMFFile(*args)

Bases: dolfin.cpp.common.Variable

Read and write Mesh , Function , MeshFunction and other objects in XDMF. This class supports the output of meshes and functions in XDMF (http://www.xdmf.org ) format. It creates an XML file that describes the data and points to a HDF5 file that stores the actual problem data. Output of data in parallel is supported. XDMF is not suitable for checkpointing as it may decimate some data.

Constructor.

Parameters:
  • comm (MPI_Comm) –
  • std::string filename (const) –
Encoding_ASCII = 1
Encoding_HDF5 = 0
close()

Close the file This closes any open HDF5 files. In ASCII mode the XML file is closed each time it is written to or read from, so close has no effect. From Python you can also use XDMFFile as a context manager:

with XDMFFile(mpi_comm_world(), 'name.xdmf') as xdmf:
    xdmf.write(mesh)

The file is automatically closed at the end of the with block

Return type:void
default_encoding = 0
read()

Read MeshValueCollection from file, optionally specifying dataset name

Parameters:
  • double > & mvc (MeshValueCollection<) – (MeshValueCollection<double>) MeshValueCollection to restore
  • name (std::string) – (std::string) Name of data attribute in XDMF file
Return type:

void

read_checkpoint()

Read a function from the XDMF file. Supplied function must come with already initialized and compatible function space. Functions saved as time-series must be addressed through (integer) internal counter and function name. Function name is a name of function that was saved. Counter stands for the number of function from time-series, e.g. counter=0 refers to first saved function regardless of its time-step value.

Parameters:
  • & u (Function) – (Function ) A function to read.
  • func_name (std::string) – (string) A name of a function to read. Must be the same on all processes in parallel.
  • counter (std::int64_t) – (int64_t) Internal integer counter - used in time-series. Default value is -1 which points to last saved function. Counter works same as python array position key, i.e. counter = -2 points to the function before the last one.
Return type:

void

thisown

The membership flag

write()

Save a cloud of points, with scalar values using an associated HDF5 file, or storing the data inline as XML.

Parameters:
  • std::vector< Point > & points (const) – (std::vector<Point>) A list of points to save.
  • std::vector< double > & values (const) – (std::vector<double>) A list of values at each point.
  • encoding (Encoding) – (Encoding) Encoding to use: HDF5 or ASCII
Return type:

void

write_checkpoint()

Save a Function to XDMF file for checkpointing, using an associated HDF5 file, or storing the data inline as XML. If the file where we would like to write exists, then the function is appended to the file.

Parameters:
  • Function & u (const) – (Function ) A function to save.
  • function_name (std::string) – (string) Name of the function used in XDMF XML and HDF file (if encoding = HDF). The string is used to fill value for XML attribute Name of Grid node. It is also used in HDF file in path to datasets. Must be the same on all processes in parallel.
  • time_step (double) – (double) Time step. It is saved only in XDMF file. The function could be saved with the same time step several times. There is an internal “counter” value stored in XDMF which differentiates the same time steps.
  • encoding (Encoding) – (Encoding) Encoding to use: HDF5 or ASCII
Return type:

void