![]() |
AlmaBTE
1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
|
Definitions corresponding to shengbte_iter.hpp. More...
#include <shengbte_iter.hpp>
Functions | |
Eigen::MatrixXd | alma::calc_shengbte_kappa (const alma::Crystal_structure &poscar, const alma::Gamma_grid &grid, const alma::Symmetry_operations &syms, std::vector< alma::Threeph_process > &threeph_procs, std::vector< alma::Twoph_process > &twoph_procs, double T, const boost::mpi::communicator &comm, double tolerance=1e-4, std::size_t maxiter=1000) |
Compute the full thermal conductivity tensor use the Omini-Sparavigna iterative approach. More... | |
Definitions corresponding to shengbte_iter.hpp.
Eigen::MatrixXd alma::calc_shengbte_kappa | ( | const alma::Crystal_structure & | poscar, |
const alma::Gamma_grid & | grid, | ||
const alma::Symmetry_operations & | syms, | ||
std::vector< alma::Threeph_process > & | threeph_procs, | ||
std::vector< alma::Twoph_process > & | twoph_procs, | ||
double | T, | ||
const boost::mpi::communicator & | comm, | ||
double | tolerance = 1e-4 , |
||
std::size_t | maxiter = 1000 |
||
) |
Compute the full thermal conductivity tensor use the Omini-Sparavigna iterative approach.
[in] | poscar | - description of the unit cell |
[in] | grid | - phonon spectrum on a regular q-point grid |
[in] | syms | - symmetry operations object |
[in,out] | threeph_procs | - list of 3-phonon processes |
[in,out] | twoph_procs | - list of 2-phonon processes |
[in] | T | - temperature in K |
[in] | comm | - communicator used to synchronize all processes |
[in] | tolerance | - maximum change in the norm betwwen iterations used as the convergence criterion |
[in] | maxiter | - maximum number of iterations before giving up |