AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
shengbte_iter.hpp
Go to the documentation of this file.
1 // Copyright 2015-2018 The ALMA Project Developers
2 //
3 // Licensed under the Apache License, Version 2.0 (the "License");
4 // you may not use this file except in compliance with the License.
5 // You may obtain a copy of the License at
6 //
7 // http://www.apache.org/licenses/LICENSE-2.0
8 //
9 // Unless required by applicable law or agreed to in writing, software
10 // distributed under the License is distributed on an "AS IS" BASIS,
11 // WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or
12 // implied. See the License for the specific language governing
13 // permissions and limitations under the License.
14 
15 #pragma once
16 
20 
21 #include <Eigen/Dense>
22 #include <processes.hpp>
23 #include <isotopic_scattering.hpp>
24 
25 namespace alma {
29 private:
31  const std::size_t nqpoints;
33  const std::size_t nirred;
35  const std::size_t nbranches;
37  const double V;
39  std::size_t n;
41  Eigen::MatrixXd omega;
43  Eigen::MatrixXd vg;
45  Eigen::MatrixXd tau0;
48  Eigen::MatrixXd F;
50  const boost::mpi::communicator& comm;
51 
52 public:
63  const Gamma_grid& grid,
64  const Symmetry_operations& syms,
65  std::vector<Threeph_process>& threeph_procs,
66  std::vector<Twoph_process>& twoph_procs,
67  double T,
68  const boost::mpi::communicator& comm_);
69 
70 
79  void next(const Crystal_structure& poscar,
80  const Gamma_grid& grid,
81  const Symmetry_operations& syms,
82  std::vector<Threeph_process>& threeph_procs,
83  std::vector<Twoph_process>& twoph_procs,
84  double T);
85 
86 
92  Eigen::Matrix3d calc_current_kappa(double T) const;
93 
94 
105  Eigen::Matrix3d calc_current_kappa_branch(double T,
106  std::size_t branch) const;
107 
108 
116  std::vector<Eigen::Matrix3d> calc_cumulative_kappa_omega(
117  double T,
118  Eigen::ArrayXd ticks);
119 
120 
128  std::vector<Eigen::Matrix3d> calc_cumulative_kappa_lambda(
129  double T,
130  Eigen::ArrayXd ticks);
131 
132 
143  Eigen::ArrayXXd calc_w() const;
144 
145 
156  Eigen::ArrayXXd calc_lambda() const;
157 };
158 
159 
175 Eigen::MatrixXd calc_shengbte_kappa(
176  const alma::Crystal_structure& poscar,
177  const alma::Gamma_grid& grid,
178  const alma::Symmetry_operations& syms,
179  std::vector<alma::Threeph_process>& threeph_procs,
180  std::vector<alma::Twoph_process>& twoph_procs,
181  double T,
182  const boost::mpi::communicator& comm,
183  double tolerance = 1e-4,
184  std::size_t maxiter = 1000);
185 } // namespace alma
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&#39;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