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

Representation of a three-phonon process. More...

#include <processes.hpp>

Collaboration diagram for alma::Threeph_process:

Public Member Functions

 Threeph_process (std::size_t _c, const std::array< std::size_t, 3 > &_q, const std::array< std::size_t, 3 > &_alpha, threeph_type _type, double _domega, double _sigma)
 Basic constructor.
 
 Threeph_process (const Threeph_process &original)
 Copy constructor.
 
double compute_gaussian ()
 Compute and return the Gaussian factor of the process, i.e., the amplitude of the regularized delta. More...
 
double compute_weighted_gaussian (const Gamma_grid &grid, double T)
 Compute and return a weighted version of the Gaussian factor of the process, containing essentially the same ingredients as the scattering rate except for the matrix element and some constants. More...
 
double compute_vp2 (const Crystal_structure &cell, const Gamma_grid &grid, const std::vector< Thirdorder_ifcs > &thirdorder)
 Compute, return and store the modulus squared of the matrix element of the process. More...
 
double get_vp2 () const
 Return the modulus squared of the matrix element of the process or throw an exception if it has not been calculated. More...
 
bool is_vp2_computed () const
 
double compute_gamma (const Gamma_grid &grid, double T)
 Get the "partial scattering rate" Gamma for this process. More...
 
double compute_gamma_reduced (const Gamma_grid &grid, double T)
 Get the reduced "partial scattering rate", meaning Gamma without the scaling factor involving the BE occupations. More...
 
Eigen::ArrayXd compute_collision (const Gamma_grid &grid, const Eigen::ArrayXXd &n0)
 Get the coefficients of the three contributions to the linearized scattering operator from this process. More...
 

Public Attributes

const std::size_t c
 Equivalence class of the first phonon.
 
const std::array< std::size_t, 3 > q
 q point indices of each of the three phonons involved.
 
const std::array< std::size_t, 3 > alpha
 Mode indices of the three phonons involved.
 
const threeph_type type
 Type of process.
 

Friends

class boost::serialization::access
 
template<class Archive >
void boost::serialization::save_construct_data (Archive &ar, const Threeph_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 three-phonon process.

Member Function Documentation

◆ compute_collision()

Eigen::ArrayXd alma::Threeph_process::compute_collision ( const Gamma_grid grid,
const Eigen::ArrayXXd &  n0 
)

Get the coefficients of the three contributions to the linearized scattering operator from this process.

vp2 must have been precomputed.

Parameters
[in]grid- regular grid with phonon spectrum
[in]n0- precomputed Bose-Einstein distribution
Returns
an array with the three coefficients

◆ compute_gamma()

double alma::Threeph_process::compute_gamma ( const Gamma_grid grid,
double  T 
)

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

vp2 must have been precomputed.

Parameters
[in]grid- regular grid with phonon spectrum
[in]T- temperature in K
Returns
Gamma, the partial scattering rate

◆ compute_gamma_reduced()

double alma::Threeph_process::compute_gamma_reduced ( const Gamma_grid grid,
double  T 
)

Get the reduced "partial scattering rate", meaning Gamma without the scaling factor involving the BE occupations.

vp2 must have been precomputed.

Parameters
[in]grid- regular grid with phonon spectrum
[in]T- temperature in K
Returns
Gamma, the partial scattering rate

◆ compute_gaussian()

double alma::Threeph_process::compute_gaussian ( )
inline

Compute and return the Gaussian factor of the process, i.e., the amplitude of the regularized delta.

The result is cached. It can be used to obtain the phase space volume of three-phonon processes.

Returns
the Gaussian factor of the process

◆ compute_vp2()

double alma::Threeph_process::compute_vp2 ( const Crystal_structure cell,
const Gamma_grid grid,
const std::vector< Thirdorder_ifcs > &  thirdorder 
)

Compute, return and store the modulus squared of the matrix element of the process.

Parameters
[in]cell- description of the unit cell
[in]grid- regular grid with phonon spectrum
[in]thirdorder- third-order ifcs
Returns
the modulus squared of the matrix element of the process

◆ compute_weighted_gaussian()

double alma::Threeph_process::compute_weighted_gaussian ( const Gamma_grid grid,
double  T 
)

Compute and return a weighted version of the Gaussian factor of the process, containing essentially the same ingredients as the scattering rate except for the matrix element and some constants.

The result can be used to obtain the weighted phase space volume of three-phonon processes. The Gaussian factor is cached.

Parameters
[in]grid- regular grid with phonon spectrum
[in]T- temperature in K
Returns
the Bose-Einstein weighted Gaussian factor of the process at the given temperature.

◆ get_vp2()

double alma::Threeph_process::get_vp2 ( ) const
inline

Return the modulus squared of the matrix element of the process or throw an exception if it has not been calculated.

Returns
the modulus squared of the matrix element of the process

◆ is_vp2_computed()

bool alma::Threeph_process::is_vp2_computed ( ) const
inline
Returns
the value of vp2_computed.

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: