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

Code implementing isotopic scattering according to Tamura's formula: S. More...

#include <cstddef>
#include <cmath>
#include <array>
#include <vector>
#include <boost/mpi.hpp>
#include <boost/math/distributions/normal.hpp>
#include <boost/serialization/array.hpp>
#include <structures.hpp>
#include <qpoint_grid.hpp>
Include dependency graph for isotopic_scattering.hpp:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

class  alma::Twoph_process
 Representation of a elastic two-phonon process. More...
 

Functions

template<class Archive >
void boost::serialization::save_construct_data (Archive &ar, const alma::Twoph_process *t, const unsigned int file_version)
 Overload required to serialize the const members of alma::Twoph_process. 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...
 
template<class Archive >
void boost::serialization::load_construct_data (Archive &ar, alma::Twoph_process *t, const unsigned int file_version)
 Overload required to deserialize the const members of alma::Twoph_process by calling the non-default constructor in place. More...
 

Detailed Description

Code implementing isotopic scattering according to Tamura's formula: S.

-I. Tamura, Isotope scattering of dispersive phonons in Ge, Phys. Rev. B 27 (1983) 858–866 A. Kundu, N. Mingo, D.A. Broido, D.A. Stewart, Role of light and heavy embedded nanoparticles on the thermal conductivity of SiGe alloys, Phys. Rev. B 84 (2011) 125426. The same method is used for alloys in the virtual crystal approximation. The classes and functions in this file are modelled after those declared in processes.hpp. There are, however, some differences. The most important aspects are:

  • Conservation of momentum is not enforced for two-phonon processes.
  • The Gaussian factor of each process is computed when the object is built.
  • The matrix element of each process is not stored in the object.

Function Documentation

◆ 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

◆ load_construct_data()

template<class Archive >
void boost::serialization::load_construct_data ( Archive &  ar,
alma::Twoph_process t,
const unsigned int  file_version 
)
inline

Overload required to deserialize the const members of alma::Twoph_process by calling the non-default constructor in place.

See the boost::serialization documentation for details.

◆ save_construct_data()

template<class Archive >
void boost::serialization::save_construct_data ( Archive &  ar,
const alma::Twoph_process t,
const unsigned int  file_version 
)
inline

Overload required to serialize the const members of alma::Twoph_process.

See the boost::serialization documentation for details.