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

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...
 

Detailed Description

Class implementing an iterative solution to the BTE following the scheme devised by Omini and Sparavigna.

Constructor & Destructor Documentation

◆ ShengBTE_iterator()

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.

Parameters
[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

Member Function Documentation

◆ calc_cumulative_kappa_lambda()

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.

Parameters
[in]T- temperature in K

◆ calc_cumulative_kappa_omega()

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.

Parameters
[in]T- temperature in K

◆ calc_current_kappa()

Eigen::Matrix3d alma::ShengBTE_iterator::calc_current_kappa ( double  T) const

Get the current estimate of the thermal conductivity tensor.

Parameters
[in]T- temperature in K
Returns
a 3x3 matrix with all the components of kappa [W / (m K)]

◆ calc_current_kappa_branch()

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.

Parameters
[in]T- temperature in K
[in]branch- branch index
Returns
a 3x3 matrix with all the components of kappa [W / (m K)]

◆ calc_lambda()

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.

Returns
an array of pseudo-scattering rates for each q point and each mode [ps ** (-1)]

◆ calc_w()

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.

Returns
an array of pseudo-scattering rates for each q point and each mode [ps ** (-1)]

◆ next()

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.

Parameters
[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

The documentation for this class was generated from the following files: