21 #include <Eigen/Dense> 31 const std::size_t nqpoints;
33 const std::size_t nirred;
35 const std::size_t nbranches;
41 Eigen::MatrixXd omega;
50 const boost::mpi::communicator& comm;
65 std::vector<Threeph_process>& threeph_procs,
66 std::vector<Twoph_process>& twoph_procs,
68 const boost::mpi::communicator& comm_);
82 std::vector<Threeph_process>& threeph_procs,
83 std::vector<Twoph_process>& twoph_procs,
106 std::size_t branch)
const;
118 Eigen::ArrayXd ticks);
130 Eigen::ArrayXd ticks);
143 Eigen::ArrayXXd
calc_w()
const;
179 std::vector<alma::Threeph_process>& threeph_procs,
180 std::vector<alma::Twoph_process>& twoph_procs,
182 const boost::mpi::communicator& comm,
183 double tolerance = 1e-4,
184 std::size_t maxiter = 1000);
Definition: analytic1d.hpp:26
Eigen::Matrix3d calc_current_kappa(double T) const
Get the current estimate of the thermal conductivity tensor.
Definition: shengbte_iter.cpp:148
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.
Definition: shengbte_iter.cpp:20
Eigen::MatrixXd 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.
Definition: shengbte_iter.cpp:266
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.
Definition: shengbte_iter.cpp:66
Class implementing an iterative solution to the BTE following the scheme devised by Omini and Sparavi...
Definition: shengbte_iter.hpp:28
Eigen::ArrayXXd calc_w() const
Obtain a set of pseudo-scattering rates from the current estimate of the solution.
Definition: shengbte_iter.cpp:179
Code implementing isotopic scattering according to Tamura's formula: S.
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-...
Definition: shengbte_iter.cpp:236
Detection and representation of allowed three-phonon processes.
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...
Definition: shengbte_iter.cpp:162
Objects of this class hold a subset of the information provided by spg_get_dataset().
Definition: symmetry.hpp:47
Objects of this class represent a regular grid with the Gamma point in one corner.
Definition: qpoint_grid.hpp:32
Hold information about a crystal structure.
Definition: structures.hpp:51
Eigen::ArrayXXd calc_lambda() const
Obtain a set of pseudo-mean free paths from the current estimate of the solution. ...
Definition: shengbte_iter.cpp:193
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...
Definition: shengbte_iter.cpp:207