AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
alma::Twoph_process Class Reference

Representation of a elastic two-phonon process. More...

#include <isotopic_scattering.hpp>

Collaboration diagram for alma::Twoph_process:

Public Member Functions

 Twoph_process (std::size_t _c, const std::array< std::size_t, 2 > &_q, const std::array< std::size_t, 2 > &_alpha, double _domega, double _sigma)
 Basic constructor.
 
 Twoph_process (const Twoph_process &original)
 Copy constructor.
 
 Twoph_process (std::size_t _c, const std::array< std::size_t, 2 > &_q, const std::array< std::size_t, 2 > &_alpha, double _domega, double _sigma, double _gaussian)
 Explicit constructor used when deserializing objects of this class. More...
 
double compute_gamma (const Crystal_structure &cell, const Gamma_grid &grid) const
 Get the "partial scattering rate" Gamma for this process. More...
 
double compute_gamma (const Crystal_structure &cell, const Eigen::Ref< const Eigen::VectorXd > &gfactors, const Gamma_grid &grid) const
 Get the "partial scattering rate" Gamma for this process. More...
 

Public Attributes

const std::size_t c
 Equivalence class of the first phonon.
 
const std::array< std::size_t, 2 > q
 q point indices of each of the two phonons involved.
 
const std::array< std::size_t, 2 > alpha
 Mode indices of the two phonons involved.
 

Friends

class boost::serialization::access
 
template<class Archive >
void boost::serialization::save_construct_data (Archive &ar, const Twoph_process *t, const unsigned int file_version)
 
void 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 > > > load_bulk_hdf5 (const char *filename, const boost::mpi::communicator &comm)
 Reconstruct all data structures from a file saved with save_bulk_hdf5(). More...
 

Detailed Description

Representation of a elastic two-phonon process.

Constructor & Destructor Documentation

◆ Twoph_process()

alma::Twoph_process::Twoph_process ( std::size_t  _c,
const std::array< std::size_t, 2 > &  _q,
const std::array< std::size_t, 2 > &  _alpha,
double  _domega,
double  _sigma,
double  _gaussian 
)
inline

Explicit constructor used when deserializing objects of this class.

Save time by avoiding the calculatio of this->gaussian.

Member Function Documentation

◆ compute_gamma() [1/2]

double alma::Twoph_process::compute_gamma ( const Crystal_structure cell,
const Gamma_grid grid 
) const

Get the "partial scattering rate" Gamma for this process.

Parameters
[in]cell- description of the unit cell
[in]grid- regular grid with phonon spectrum
Returns
Gamma, the partial scattering rate

◆ compute_gamma() [2/2]

double alma::Twoph_process::compute_gamma ( const Crystal_structure cell,
const Eigen::Ref< const Eigen::VectorXd > &  gfactors,
const Gamma_grid grid 
) const

Get the "partial scattering rate" Gamma for this process.

using a custom set of g factors.

Parameters
[in]cell- description of the unit cell
[in]gfactors- Pearson deviation coefficient of the mass at each site.
[in]grid- regular grid with phonon spectrum
Returns
Gamma, the partial scattering rate

Friends And Related Function Documentation

◆ 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> > > load_bulk_hdf5 ( const char *  filename,
const boost::mpi::communicator &  comm 
)
friend

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

◆ save_bulk_hdf5

void 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 
)
friend

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

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