AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
green1d.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 <iostream>
22 #include <fstream>
23 #include <string>
24 #include <vector>
25 #include <boost/filesystem.hpp>
26 #include <dynamical_matrix.hpp>
27 #include <qpoint_grid.hpp>
28 #include <vasp_io.hpp>
29 
30 namespace alma {
35 public:
37  const std::size_t nqpoints;
39  const std::size_t ndof;
50  Green1d_factory(const Crystal_structure& structure,
51  const Eigen::Ref<const Eigen::VectorXd>& q0,
52  const Eigen::Ref<const Eigen::VectorXi>& normal,
53  const std::size_t nqpoints_,
54  const Dynamical_matrix_builder& builder);
61  Eigen::MatrixXcd build(double omega, std::size_t ncells) const;
62 
73  double calc_w0_dmass(
74  double q,
75  double omega,
76  const Eigen::Ref<const Eigen::VectorXcd>& wfin,
77  const Eigen::Ref<const Eigen::VectorXd>& factors) const;
78 
79 private:
82  const double offset;
85  Eigen::ArrayXd qgrid;
89  Eigen::ArrayXXd Egrid;
94  Eigen::ArrayXXd dEgrid;
98  std::vector<Eigen::MatrixXcd> wfgrid;
105  inline Eigen::MatrixXcd get_superwfs(std::size_t iq,
106  std::size_t ncells) const;
107 
113  Eigen::MatrixXcd compute_weights(double omega2) const;
114 };
115 } // namespace alma
const std::size_t nqpoints
Number of q points in the grid.
Definition: green1d.hpp:37
Definition: analytic1d.hpp:26
Eigen::MatrixXcd build(double omega, std::size_t ncells) const
Obtain the Green&#39;s function for a particular angular frequency.
Definition: green1d.cpp:57
Objects of this class enable the calling code to compute 1d Green&#39;s functions along particular direct...
Definition: green1d.hpp:34
Code related to the dynamical matrix.
Classes and functions used to manipulate grids in reciprocal space.
const std::size_t ndof
Number of degrees of freedom in the unit cell.
Definition: green1d.hpp:39
Factory of Dynamical_matrix objects.
Definition: dynamical_matrix.hpp:108
Green1d_factory(const Crystal_structure &structure, const Eigen::Ref< const Eigen::VectorXd > &q0, const Eigen::Ref< const Eigen::VectorXi > &normal, const std::size_t nqpoints_, const Dynamical_matrix_builder &builder)
Basic constructor.
Definition: green1d.cpp:30
Routines to load files coming from the VASP + Phonopy ecosystem.
Hold information about a crystal structure.
Definition: structures.hpp:51
double calc_w0_dmass(double q, double omega, const Eigen::Ref< const Eigen::VectorXcd > &wfin, const Eigen::Ref< const Eigen::VectorXd > &factors) const
Obtain the scattering rate caused by a diagonal perturbation proportional to the frequency squared...
Definition: green1d.cpp:76