DOLFIN
DOLFIN C++ interface
Public Types | Public Member Functions | Static Public Attributes | List of all members
dolfin::XDMFFile Class Reference

Read and write Mesh, Function, MeshFunction and other objects in XDMF. More...

#include <XDMFFile.h>

Inheritance diagram for dolfin::XDMFFile:
Inheritance graph
[legend]
Collaboration diagram for dolfin::XDMFFile:
Collaboration graph
[legend]

Public Types

enum  Encoding { HDF5, ASCII }
 File encoding type.
 

Public Member Functions

 XDMFFile (const std::string filename)
 Constructor.
 
 XDMFFile (MPI_Comm comm, const std::string filename)
 Constructor.
 
 ~XDMFFile ()
 Destructor.
 
void close ()
 
void write (const Mesh &mesh, Encoding encoding=default_encoding)
 
void write_checkpoint (const Function &u, std::string function_name, double time_step=0.0, Encoding encoding=default_encoding, bool append=false)
 
void write (const Function &u, Encoding encoding=default_encoding)
 
void write (const Function &u, double t, Encoding encoding=default_encoding)
 
void write (const MeshFunction< bool > &meshfunction, Encoding encoding=default_encoding)
 
void write (const MeshFunction< int > &meshfunction, Encoding encoding=default_encoding)
 
void write (const MeshFunction< std::size_t > &meshfunction, Encoding encoding=default_encoding)
 
void write (const MeshFunction< double > &meshfunction, Encoding encoding=default_encoding)
 
void write (const MeshValueCollection< bool > &mvc, Encoding encoding=default_encoding)
 
void write (const MeshValueCollection< int > &mvc, Encoding encoding=default_encoding)
 
void write (const MeshValueCollection< std::size_t > &mvc, Encoding encoding=default_encoding)
 
void write (const MeshValueCollection< double > &mvc, Encoding encoding=default_encoding)
 
void write (const std::vector< Point > &points, Encoding encoding=default_encoding)
 
void write (const std::vector< Point > &points, const std::vector< double > &values, Encoding encoding=default_encoding)
 
void read (Mesh &mesh) const
 
void read_checkpoint (Function &u, std::string func_name, std::int64_t counter=-1)
 
void read (MeshFunction< bool > &meshfunction, std::string name="")
 
void read (MeshFunction< int > &meshfunction, std::string name="")
 
void read (MeshFunction< std::size_t > &meshfunction, std::string name="")
 
void read (MeshFunction< double > &meshfunction, std::string name="")
 
void read (MeshValueCollection< bool > &mvc, std::string name="")
 
void read (MeshValueCollection< int > &mvc, std::string name="")
 
void read (MeshValueCollection< std::size_t > &mvc, std::string name="")
 
void read (MeshValueCollection< double > &mvc, std::string name="")
 
template<>
void add_data_item (MPI_Comm comm, pugi::xml_node &xml_node, hid_t h5_id, const std::string h5_path, const std::vector< bool > &x, const std::vector< std::int64_t > shape, const std::string number_type)
 
- Public Member Functions inherited from dolfin::Variable
 Variable ()
 Create unnamed variable.
 
 Variable (const std::string name, const std::string label)
 Create variable with given name and label.
 
 Variable (const Variable &variable)
 Copy constructor.
 
virtual ~Variable ()
 Destructor.
 
const Variableoperator= (const Variable &variable)
 Assignment operator.
 
void rename (const std::string name, const std::string label)
 Rename variable.
 
std::string name () const
 Return name.
 
std::string label () const
 Return label (description)
 
std::size_t id () const
 
virtual std::string str (bool verbose) const
 Return informal string representation (pretty-print)
 

Static Public Attributes

static const Encoding default_encoding = Encoding::HDF5
 Default encoding type.
 

Additional Inherited Members

- Public Attributes inherited from dolfin::Variable
Parameters parameters
 Parameters.
 

Detailed Description

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.

Member Function Documentation

◆ close()

void XDMFFile::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

◆ read() [1/9]

void XDMFFile::read ( Mesh mesh) const

Read in the first Mesh in XDMF file

Parameters
mesh(Mesh) Mesh to fill from XDMF file

◆ read() [2/9]

void XDMFFile::read ( MeshFunction< bool > &  meshfunction,
std::string  name = "" 
)

Read first MeshFunction from file

Parameters
meshfunction(MeshFunction<bool>) MeshFunction to restore
name(std::string) Name of data attribute in XDMF file

◆ read() [3/9]

void XDMFFile::read ( MeshFunction< int > &  meshfunction,
std::string  name = "" 
)

Read first MeshFunction from file

Parameters
meshfunction(MeshFunction<int>) MeshFunction to restore
name(std::string) Name of data attribute in XDMF file

◆ read() [4/9]

void XDMFFile::read ( MeshFunction< std::size_t > &  meshfunction,
std::string  name = "" 
)

Read MeshFunction from file, optionally specifying dataset name

Parameters
meshfunction(MeshFunction<std::size_t>) MeshFunction to restore
name(std::string) Name of data attribute in XDMF file

◆ read() [5/9]

void XDMFFile::read ( MeshFunction< double > &  meshfunction,
std::string  name = "" 
)

Read MeshFunction from file, optionally specifying dataset name

Parameters
meshfunction(MeshFunction<double>) MeshFunction to restore
name(std::string) Name of data attribute in XDMF file

◆ read() [6/9]

void XDMFFile::read ( MeshValueCollection< bool > &  mvc,
std::string  name = "" 
)

Read MeshValueCollection from file, optionally specifying dataset name

Parameters
mvc(MeshValueCollection<bool>) MeshValueCollection to restore
name(std::string) Name of data attribute in XDMF file

◆ read() [7/9]

void XDMFFile::read ( MeshValueCollection< int > &  mvc,
std::string  name = "" 
)

Read MeshValueCollection from file, optionally specifying dataset name

Parameters
mvc(MeshValueCollection<int>) MeshValueCollection to restore
name(std::string) Name of data attribute in XDMF file

◆ read() [8/9]

void XDMFFile::read ( MeshValueCollection< std::size_t > &  mvc,
std::string  name = "" 
)

Read MeshValueCollection from file, optionally specifying dataset name

Parameters
mvc(MeshValueCollection<std::size_t>) MeshValueCollection to restore
name(std::string) Name of data attribute in XDMF file

◆ read() [9/9]

void XDMFFile::read ( MeshValueCollection< double > &  mvc,
std::string  name = "" 
)

Read MeshValueCollection from file, optionally specifying dataset name

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

◆ read_checkpoint()

void XDMFFile::read_checkpoint ( Function u,
std::string  func_name,
std::int64_t  counter = -1 
)

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) A function to read.
func_name(string) A name of a function to read. Must be the same on all processes in parallel.
counter(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.

◆ write() [1/13]

void XDMFFile::write ( const Mesh mesh,
Encoding  encoding = default_encoding 
)

Save a mesh to XDMF format, either using an associated HDF5 file, or storing the data inline as XML Create function on given function space

Parameters
mesh(Mesh) A mesh to save.
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [2/13]

void XDMFFile::write ( const Function u,
Encoding  encoding = default_encoding 
)

Save a Function to XDMF file for visualisation, using an associated HDF5 file, or storing the data inline as XML.

Parameters
u(Function) A function to save.
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [3/13]

void XDMFFile::write ( const Function u,
double  t,
Encoding  encoding = default_encoding 
)

Save a Function with timestamp to XDMF file for visualisation, using an associated HDF5 file, or storing the data inline as XML.

You can control the output with the following boolean parameters on the XDMFFile class:

  • rewrite_function_mesh (default true): Controls whether the mesh will be rewritten every timestep. If the mesh does not change this can be turned off to create smaller files.
  • functions_share_mesh (default false): Controls whether all functions on a single time step share the same mesh. If true the files created will be smaller and also behave better in Paraview, at least in version 5.3.0
Parameters
u(Function) A function to save.
t(double) Timestep
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [4/13]

void XDMFFile::write ( const MeshFunction< bool > &  meshfunction,
Encoding  encoding = default_encoding 
)

Save MeshFunction to file using an associated HDF5 file, or storing the data inline as XML.

Parameters
meshfunction(MeshFunction) A meshfunction to save.
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [5/13]

void XDMFFile::write ( const MeshFunction< int > &  meshfunction,
Encoding  encoding = default_encoding 
)

Save MeshFunction to file using an associated HDF5 file, or storing the data inline as XML.

Parameters
meshfunction(MeshFunction) A meshfunction to save.
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [6/13]

void XDMFFile::write ( const MeshFunction< std::size_t > &  meshfunction,
Encoding  encoding = default_encoding 
)

Save MeshFunction to file using an associated HDF5 file, or storing the data inline as XML.

Parameters
meshfunction(MeshFunction) A meshfunction to save.
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [7/13]

void XDMFFile::write ( const MeshFunction< double > &  meshfunction,
Encoding  encoding = default_encoding 
)

Save MeshFunction to file using an associated HDF5 file, or storing the data inline as XML.

Parameters
meshfunction(MeshFunction) A meshfunction to save.
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [8/13]

void XDMFFile::write ( const MeshValueCollection< bool > &  mvc,
Encoding  encoding = default_encoding 
)

Write out mesh value collection (subset) using an associated HDF5 file, or storing the data inline as XML.

Parameters
mvc(MeshValueCollection<bool>) MeshValueCollection to save
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [9/13]

void XDMFFile::write ( const MeshValueCollection< int > &  mvc,
Encoding  encoding = default_encoding 
)

Write out mesh value collection (subset) using an associated HDF5 file, or storing the data inline as XML.

Parameters
mvc(MeshValueCollection<int>) MeshValueCollection to save
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [10/13]

void XDMFFile::write ( const MeshValueCollection< std::size_t > &  mvc,
Encoding  encoding = default_encoding 
)

Write out mesh value collection (subset) using an associated HDF5 file, or storing the data inline as XML.

Parameters
mvc(MeshValueCollection<int>) MeshValueCollection to save
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [11/13]

void XDMFFile::write ( const MeshValueCollection< double > &  mvc,
Encoding  encoding = default_encoding 
)

Write out mesh value collection (subset) using an associated HDF5 file, or storing the data inline as XML.

Parameters
mvc(MeshValueCollection<double>) MeshValueCollection to save
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [12/13]

void XDMFFile::write ( const std::vector< Point > &  points,
Encoding  encoding = default_encoding 
)

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

Parameters
points(std::vector<Point>) A list of points to save.
encoding(Encoding) Encoding to use: HDF5 or ASCII

◆ write() [13/13]

void XDMFFile::write ( const std::vector< Point > &  points,
const std::vector< double > &  values,
Encoding  encoding = default_encoding 
)

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

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

◆ write_checkpoint()

void XDMFFile::write_checkpoint ( const Function u,
std::string  function_name,
double  time_step = 0.0,
Encoding  encoding = default_encoding,
bool  append = false 
)

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
u(Function) A function to save.
function_name(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) 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 to use: HDF5 or ASCII
append(default false) If set (default) to false, then XDMF and HDF files (if encoding = HDF) are overwritten. If set to true, then both XDMF and HDF (if encoding = HDF) files are appended at end. Set append to true when saving timeseries.

The documentation for this class was generated from the following files: