![]() |
AlmaBTE
1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
|
Definitions corresponding to bulk_hdf5.hpp. More...

Functions | |
| void | alma::write_string_attribute (H5Location &loc, const std::string &name, const std::string &value) |
| Add a string-valued attribute to an HDF5 file, group, or dataset. More... | |
| template<typename T > | |
| void | alma::create_eigen_dataset (CommonFG &loc, const Eigen::DenseBase< T > &data, const std::string &name, const std::string &units) |
| Create a dataset from an Eigen object. More... | |
| void | alma::write_file_attributes (H5File &file, const std::string &description) |
| Write some metadata to the HDF5 file as attributes. More... | |
| void | alma::create_groups (H5File &file) |
| Create some predefined groups in the HDF5 file. More... | |
| void | alma::save_bulk_hdf5 (const char *filename, const std::string &description, const Crystal_structure &cell, const Symmetry_operations &symmetries, const Gamma_grid &grid, const std::vector< Threeph_process > &processes, const boost::mpi::communicator &comm) |
| Save all relevant information about a bulk material to an HDF5 file. More... | |
| std::string | alma::read_file_attributes (H5File &file) |
| Read an HDF5 file's metadata and check if the versions are compatible. More... | |
| std::tuple< std::string, std::unique_ptr< Crystal_structure >, std::unique_ptr< Symmetry_operations >, std::unique_ptr< Gamma_grid >, std::unique_ptr< std::vector< Threeph_process > > > | alma::load_bulk_hdf5 (const char *filename, const boost::mpi::communicator &comm) |
| Reconstruct all data structures from a file saved with save_bulk_hdf5(). More... | |
| std::vector< std::string > | alma::list_scattering_subgroups (const char *filename, const boost::mpi::communicator &comm) |
| Return the names of all subgroups in the "/scattering" group of an HDF5 file. More... | |
| Scattering_subgroup | alma::load_scattering_subgroup (const char *filename, const std::string &groupname, const boost::mpi::communicator &comm) |
| Read the scattering rates stored in a subgroup from the "/scattering" group. More... | |
| void | alma::write_scattering_subgroup (const char *filename, const Scattering_subgroup &group, const boost::mpi::communicator &comm) |
| Add a subgroup to the "/scattering" group of the HDF5 file. More... | |
Variables | |
| constexpr std::size_t | ABSORPTION_BIT = 0 |
| Constants used to build the "/threeph_processes/type" dataset. More... | |
| constexpr std::size_t | GAUSSIAN_BIT = 1 |
| constexpr std::size_t | VP2_BIT = 2 |
Definitions corresponding to bulk_hdf5.hpp.
| void alma::create_eigen_dataset | ( | CommonFG & | loc, |
| const Eigen::DenseBase< T > & | data, | ||
| const std::string & | name, | ||
| const std::string & | units | ||
| ) |
Create a dataset from an Eigen object.
All objects are assumed to be 2D, and their contents are assumed to be double-precission floating point numbers.
| [in,out] | loc | - the HDF5 file or group |
| [in] | data | - the Eigen object |
| [in] | name | - the name of the dataset |
| [in] | units | - the units of the data |
| void alma::create_groups | ( | H5File & | file | ) |
Create some predefined groups in the HDF5 file.
| [in,out] | file | - the hdf5 file |
| std::vector< std::string > alma::list_scattering_subgroups | ( | const char * | filename, |
| const boost::mpi::communicator & | comm | ||
| ) |
Return the names of all subgroups in the "/scattering" group of an HDF5 file.
An empty vector is returned both when the "/scattering" group does not exist and when it exists but does not contain any subset.
| [in] | filename | - path to the HDF5 file |
| [in] | comm | - MPI communicator used to coordinate with all other processes |
| std::tuple< std::string, std::unique_ptr< Crystal_structure >, std::unique_ptr< Symmetry_operations >, std::unique_ptr< Gamma_grid >, std::unique_ptr< std::vector< Threeph_process > > > alma::load_bulk_hdf5 | ( | const char * | filename, |
| const boost::mpi::communicator & | comm | ||
| ) |
Reconstruct all data structures from a file saved with save_bulk_hdf5().
| [in] | filename | - path to the HDF5 file |
| [in] | comm | - MPI communicator used to determine which part of the data to load in this process. |
| Scattering_subgroup alma::load_scattering_subgroup | ( | const char * | filename, |
| const std::string & | groupname, | ||
| const boost::mpi::communicator & | comm | ||
| ) |
Read the scattering rates stored in a subgroup from the "/scattering" group.
The existence of the group is not checked.
| [in] | filename | - path to the HDF5 file |
| [in] | groupname | - name of the group |
| [in] | comm | - MPI communicator used to coordinate with all other processes |
| std::string alma::read_file_attributes | ( | H5File & | file | ) |
Read an HDF5 file's metadata and check if the versions are compatible.
| [in] | file | - the HDF5 file |
| void alma::save_bulk_hdf5 | ( | const char * | filename, |
| const std::string & | description, | ||
| const Crystal_structure & | cell, | ||
| const Symmetry_operations & | symmetries, | ||
| const Gamma_grid & | grid, | ||
| const std::vector< Threeph_process > & | processes, | ||
| const boost::mpi::communicator & | comm | ||
| ) |
Save all relevant information about a bulk material to an HDF5 file.
| [in] | filename | - path to the HDF5 file. It can be NULL for all processes except the one with rank == 0. |
| [in] | description | - a free-format description |
| [in] | cell | - description of the unit cell |
| [in] | symmetries | - symmetry operations of the unit cell |
| [in] | grid | - phonon spectrum on a refular grid |
| [in] | processes | - allowed three-phonon processes |
| [in] | comm | - MPI communicator used to coordinate with all other processes |
| void alma::write_file_attributes | ( | H5File & | file, |
| const std::string & | description | ||
| ) |
Write some metadata to the HDF5 file as attributes.
| [in,out] | file | - the HDF5 file |
| [in] | description | - a free-format description |
| void alma::write_scattering_subgroup | ( | const char * | filename, |
| const Scattering_subgroup & | group, | ||
| const boost::mpi::communicator & | comm | ||
| ) |
Add a subgroup to the "/scattering" group of the HDF5 file.
The existence of the group is not checked.
| [in] | filename | - path to the HDF5 file |
| [in] | group | - information about the group |
| [in] | comm | - MPI communicator used to coordinate with all other processes |
| void alma::write_string_attribute | ( | H5Location & | loc, |
| const std::string & | name, | ||
| const std::string & | value | ||
| ) |
Add a string-valued attribute to an HDF5 file, group, or dataset.
Note that UTF-8 encoding is assumed.
| [in,out] | loc | - HDF5 file, group, or dataset |
| [in] | name | - attribute identifier |
| [in] | value | - attribute value |
| constexpr std::size_t ABSORPTION_BIT = 0 |
Constants used to build the "/threeph_processes/type" dataset.
See below for details.