![]() |
AlmaBTE
1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
|
Class implementing an iterative solution to the BTE following the scheme devised by Omini and Sparavigna. More...
#include <shengbte_iter.hpp>
Public Member Functions | |
ShengBTE_iterator (const Crystal_structure &poscar, const Gamma_grid &grid, const Symmetry_operations &syms, std::vector< Threeph_process > &threeph_procs, std::vector< Twoph_process > &twoph_procs, double T, const boost::mpi::communicator &comm_) | |
Basic constructor to initialize F. More... | |
void | next (const Crystal_structure &poscar, const Gamma_grid &grid, const Symmetry_operations &syms, std::vector< Threeph_process > &threeph_procs, std::vector< Twoph_process > &twoph_procs, double T) |
Advance to the next iteration. More... | |
Eigen::Matrix3d | calc_current_kappa (double T) const |
Get the current estimate of the thermal conductivity tensor. More... | |
Eigen::Matrix3d | calc_current_kappa_branch (double T, std::size_t branch) const |
Get the current estimate of the contribution to the thermal conductivity tensor from a single branch. More... | |
std::vector< Eigen::Matrix3d > | calc_cumulative_kappa_omega (double T, Eigen::ArrayXd ticks) |
Obtain the cumulative histogram of contributions to the thermal conductivity as a function of angular frequency. More... | |
std::vector< Eigen::Matrix3d > | calc_cumulative_kappa_lambda (double T, Eigen::ArrayXd ticks) |
Obtain the cumulative histogram of contributions to the thermal conductivity as a function of pseudo-mean free path. More... | |
Eigen::ArrayXXd | calc_w () const |
Obtain a set of pseudo-scattering rates from the current estimate of the solution. More... | |
Eigen::ArrayXXd | calc_lambda () const |
Obtain a set of pseudo-mean free paths from the current estimate of the solution. More... | |
Class implementing an iterative solution to the BTE following the scheme devised by Omini and Sparavigna.
alma::ShengBTE_iterator::ShengBTE_iterator | ( | const Crystal_structure & | poscar, |
const Gamma_grid & | grid, | ||
const Symmetry_operations & | syms, | ||
std::vector< Threeph_process > & | threeph_procs, | ||
std::vector< Twoph_process > & | twoph_procs, | ||
double | T, | ||
const boost::mpi::communicator & | comm_ | ||
) |
Basic constructor to initialize F.
[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 |
std::vector< Eigen::Matrix3d > alma::ShengBTE_iterator::calc_cumulative_kappa_lambda | ( | double | T, |
Eigen::ArrayXd | ticks | ||
) |
Obtain the cumulative histogram of contributions to the thermal conductivity as a function of pseudo-mean free path.
[in] | T | - temperature in K |
std::vector< Eigen::Matrix3d > alma::ShengBTE_iterator::calc_cumulative_kappa_omega | ( | double | T, |
Eigen::ArrayXd | ticks | ||
) |
Obtain the cumulative histogram of contributions to the thermal conductivity as a function of angular frequency.
[in] | T | - temperature in K |
Eigen::Matrix3d alma::ShengBTE_iterator::calc_current_kappa | ( | double | T | ) | const |
Get the current estimate of the thermal conductivity tensor.
[in] | T | - temperature in K |
Eigen::Matrix3d alma::ShengBTE_iterator::calc_current_kappa_branch | ( | double | T, |
std::size_t | branch | ||
) | const |
Get the current estimate of the contribution to the thermal conductivity tensor from a single branch.
In this context, a branch is defined as the set of modes with the same index when energies are sorted in ascending order at each q point.
[in] | T | - temperature in K |
[in] | branch | - branch index |
Eigen::ArrayXXd alma::ShengBTE_iterator::calc_lambda | ( | ) | const |
Obtain a set of pseudo-mean free paths from the current estimate of the solution.
Note that these are not real MFPs, since no such thing exist in the full linearized BTE formalism for phonons. Substituting these into the RTA expression for kappa will not yield the right thermal conductivity. Moreover, some elements can be negative.See the ShengBTE paper for details.
Eigen::ArrayXXd alma::ShengBTE_iterator::calc_w | ( | ) | const |
Obtain a set of pseudo-scattering rates from the current estimate of the solution.
Note that these are not the inverse of any real relaxation time, since no such thing exist in the full linearized BTE formalism for phonons. Substituting these into the RTA expression for kappa will not yield the right thermal conductivity. Moreover, some elements can be negative. See the ShengBTE paper for details.
void alma::ShengBTE_iterator::next | ( | const Crystal_structure & | poscar, |
const Gamma_grid & | grid, | ||
const Symmetry_operations & | syms, | ||
std::vector< Threeph_process > & | threeph_procs, | ||
std::vector< Twoph_process > & | twoph_procs, | ||
double | T | ||
) |
Advance to the next iteration.
[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 |