20 #ifndef __DOLFIN_HDF5FILE_H 21 #define __DOLFIN_HDF5FILE_H 29 #include <dolfin/common/MPI.h> 30 #include <dolfin/common/Variable.h> 31 #include <dolfin/geometry/Point.h> 32 #include "HDF5Attribute.h" 33 #include "HDF5Interface.h" 43 template<
typename T>
class MeshFunction;
44 template<
typename T>
class MeshValueCollection;
54 HDF5File(MPI_Comm comm,
const std::string filename,
55 const std::string file_mode);
67 void write(
const std::vector<Point>& points,
const std::string
name);
70 void write(
const std::vector<double>& values,
const std::string name);
78 const bool use_partition_from_file)
const;
81 void write(
const Mesh& mesh,
const std::string name);
85 void write(
const Mesh& mesh,
const std::size_t cell_dim,
86 const std::string name);
92 void write(
const Function& u,
const std::string name,
double timestamp);
107 void read(
Mesh& mesh,
const std::string data_path,
108 bool use_partition_from_file)
const;
120 const std::string topology_path,
121 const std::string geometry_path,
122 const int gdim ,
const CellType& cell_type,
123 const std::int64_t expected_num_global_cells,
124 const std::int64_t expected_num_global_points,
125 bool use_partition_from_file)
const;
129 const std::string name);
136 const std::string name);
143 const std::string name)
const;
150 const std::string name)
const;
157 const std::string name);
161 const std::string name);
165 const std::string name);
169 const std::string name)
const;
173 const std::string name)
const;
177 const std::string name)
const;
180 bool has_dataset(
const std::string dataset_name)
const;
192 {
return _hdf5_file_id; }
201 template <
typename T>
203 const std::string name);
206 template <
typename T>
208 const std::string name)
const;
211 template <
typename T>
212 void write_mesh_value_collection_old(
214 const std::string name);
218 template <
typename T>
220 const std::string name);
223 template <
typename T>
225 const std::string name)
const;
228 template <
typename T>
230 const std::string name)
const;
234 template <
typename T>
235 void write_data(
const std::string dataset_name,
236 const std::vector<T>& data,
237 const std::vector<std::int64_t> global_size,
249 template <
typename T>
250 void HDF5File::write_data(
const std::string dataset_name,
251 const std::vector<T>& data,
252 const std::vector<std::int64_t> global_size,
255 dolfin_assert(_hdf5_file_id > 0);
256 dolfin_assert(global_size.size() > 0);
259 std::size_t num_local_items = 1;
260 for (std::size_t i = 1; i < global_size.size(); ++i)
261 num_local_items *= global_size[i];
262 num_local_items = data.size()/num_local_items;
267 std::pair<std::size_t, std::size_t> range(offset,
268 offset + num_local_items);
273 std::string dset_name(dataset_name);
274 if (dset_name[0] !=
'/')
275 dset_name =
"/" + dataset_name;
278 range, global_size, use_mpi_io, chunking);
std::string name() const
Return name.
Definition: Variable.cpp:71
Common base class for DOLFIN variables.
Definition: Variable.h:35
void read(GenericVector &x, const std::string dataset_name, const bool use_partition_from_file) const
Definition: HDF5File.cpp:185
void set_mpi_atomicity(bool atomic)
Set the MPI atomicity.
Definition: HDF5File.cpp:1967
void close()
Close file.
Definition: HDF5File.cpp:109
bool has_dataset(const std::string dataset_name) const
Check if dataset exists in HDF5 file.
Definition: HDF5File.cpp:1948
static std::size_t global_offset(MPI_Comm comm, std::size_t range, bool exclusive)
Definition: MPI.cpp:185
~HDF5File()
Destructor.
Definition: HDF5File.cpp:104
Read and write Mesh, Function, MeshFunction and other objects in XDMF.
Definition: XDMFFile.h:77
Definition: HDF5Attribute.h:38
Definition: CellType.h:46
MPI_Comm comm() const
Return the underlying MPI_Comm object.
Definition: MPI.cpp:117
Parameters parameters
Parameters.
Definition: Variable.h:74
Definition: Function.h:65
Definition: GenericFile.h:38
void flush()
Flush buffered I/O to disk.
Definition: HDF5File.cpp:117
bool get_mpi_atomicity() const
Get the MPI atomicity.
Definition: HDF5File.cpp:1973
Definition: HDF5File.h:47
Definition: TimeSeries.h:46
HDF5File(MPI_Comm comm, const std::string filename, const std::string file_mode)
Definition: HDF5File.cpp:58
static void write_dataset(const hid_t file_handle, const std::string dataset_path, const std::vector< T > &data, const std::pair< std::int64_t, std::int64_t > range, const std::vector< std::int64_t > global_size, bool use_mpio, bool use_chunking)
This class defines a common interface for vectors.
Definition: GenericVector.h:47
void write(const std::vector< Point > &points, const std::string name)
Write points to file.
Definition: HDF5File.cpp:123