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

Definitions corresponding to isotopic_scattering.hpp. More...

#include <cmath>
#include <complex>
#include <limits>
#include <constants.hpp>
#include <periodic_table.hpp>
#include <utilities.hpp>
#include <isotopic_scattering.hpp>
Include dependency graph for isotopic_scattering.cpp:

Functions

std::array< double, 2 > alma::calc_percentiles_log (const Eigen::Ref< const Eigen::ArrayXXd > &sigma)
 Compute the 25th and 75th percentiles of log(sigma) More...
 
std::vector< Twoph_process > alma::find_allowed_twoph (const Gamma_grid &grid, const boost::mpi::communicator &communicator, double scalebroad=1.0)
 Look for allowed two-phonon processes in a regular grid. More...
 
Eigen::ArrayXXd alma::calc_w0_twoph (const Crystal_structure &cell, const alma::Gamma_grid &grid, const std::vector< alma::Twoph_process > &processes, const boost::mpi::communicator &comm)
 Compute and the two-phonon contribution to the RTA scattering rates for all vibrational modes on a grid. More...
 
Eigen::ArrayXXd alma::calc_w0_twoph (const Crystal_structure &cell, const Eigen::Ref< const Eigen::VectorXd > &gfactors, const alma::Gamma_grid &grid, const std::vector< alma::Twoph_process > &processes, const boost::mpi::communicator &comm)
 Compute and the two-phonon contribution to the RTA scattering rates for all vibrational modes on a grid, with a custom set of g factors. More...
 

Detailed Description

Definitions corresponding to isotopic_scattering.hpp.

Function Documentation

◆ calc_percentiles_log()

std::array<double, 2> alma::calc_percentiles_log ( const Eigen::Ref< const Eigen::ArrayXXd > &  sigma)

Compute the 25th and 75th percentiles of log(sigma)

Parameters
[in]sigmaan Eigen array with all broadening parameters
Returns
an array containing the 25th and 75th percentiles

◆ calc_w0_twoph() [1/2]

Eigen::ArrayXXd alma::calc_w0_twoph ( const Crystal_structure cell,
const alma::Gamma_grid grid,
const std::vector< alma::Twoph_process > &  processes,
const boost::mpi::communicator &  comm 
)

Compute and the two-phonon contribution to the RTA scattering rates for all vibrational modes on a grid.

Parameters
[in]cell- a description of the crystal structure
[in]grid- phonon spectrum on a regular grid
[in]processes- a vector of allowed two-phonon processes
[in]comm- an mpi communicator
Returns
a set of scattering rates

◆ calc_w0_twoph() [2/2]

Eigen::ArrayXXd alma::calc_w0_twoph ( const Crystal_structure cell,
const Eigen::Ref< const Eigen::VectorXd > &  gfactors,
const alma::Gamma_grid grid,
const std::vector< alma::Twoph_process > &  processes,
const boost::mpi::communicator &  comm 
)

Compute and the two-phonon contribution to the RTA scattering rates for all vibrational modes on a grid, with a custom set of g factors.

Parameters
[in]cell- a description of the crystal structure
[in]grid- phonon spectrum on a regular grid
[in]gfactors- Pearson deviation coefficient of the mass at each site.
[in]processes- a vector of allowed two-phonon processes
[in]comm- an mpi communicator
Returns
a set of scattering rates

◆ find_allowed_twoph()

std::vector< Twoph_process > alma::find_allowed_twoph ( const Gamma_grid grid,
const boost::mpi::communicator &  communicator,
double  scalebroad = 1.0 
)

Look for allowed two-phonon processes in a regular grid.

Iterate over part of the irreducible q points in the grid (trying to evenly split the equivalence classes over processes) and look for allowed two-phonon processes involving one phonon from that part and another phonon from anywhere in the grid.

Parameters
[in]grid- a regular grid containing Gamma
[in]communicator- MPI communicator to use
[in]scalebroad- factor modulating all the broadenings
Returns
a vector of Twoph_process objects