AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
bulk_hdf5.hpp File Reference

Functions to save and load data about a bulk system to an HDF5 file. More...

#include <string>
#include <vector>
#include <tuple>
#include <boost/mpi.hpp>
#include <structures.hpp>
#include <symmetry.hpp>
#include <qpoint_grid.hpp>
#include <processes.hpp>
Include dependency graph for bulk_hdf5.hpp:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

class  alma::Scattering_subgroup
 POD class describing a subgroup of the "/scattering" group of an HDF5 file. More...
 

Functions

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::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...
 

Detailed Description

Functions to save and load data about a bulk system to an HDF5 file.

Function Documentation

◆ list_scattering_subgroups()

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.

Parameters
[in]filename- path to the HDF5 file
[in]comm- MPI communicator used to coordinate with all other processes
Returns
a vector of subgroup names

◆ load_bulk_hdf5()

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().

Parameters
[in]filename- path to the HDF5 file
[in]comm- MPI communicator used to determine which part of the data to load in this process.
Returns
a tuple with all the information in the file, namely: description - a free-format description cell - description of the unit cell symmetries - symmetry operations of the unit cell grid - phonon spectrum on a refular grid processes - allowed three-phonon processes

◆ load_scattering_subgroup()

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.

Parameters
[in]filename- path to the HDF5 file
[in]groupname- name of the group
[in]comm- MPI communicator used to coordinate with all other processes
Returns
a Scattering_subgroup structure

◆ save_bulk_hdf5()

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.

Parameters
[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

◆ write_scattering_subgroup()

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.

Parameters
[in]filename- path to the HDF5 file
[in]group- information about the group
[in]comm- MPI communicator used to coordinate with all other processes